module_ra_gsfcsw.F
References to this file elsewhere.
1 !Comment the following out to turn off aerosol-radiation
2 !feedback between MOSAIC and GSFCSW. wig, 21-Feb-2005
3
4 MODULE module_ra_gsfcsw
5
6 REAL, PARAMETER, PRIVATE :: thresh=1.e-9
7 REAL, SAVE :: center_lat
8
9 ! Assign co2 and trace gases amount (units are parts/part by volumn)
10
11 REAL, PARAMETER, PRIVATE :: co2 = 300.e-6
12
13 CONTAINS
14
15 !------------------------------------------------------------------
16 ! urban related variable are added to arguments of gsfcswrad
17 !------------------------------------------------------------------
18 SUBROUTINE GSFCSWRAD(rthraten,gsw,xlat,xlong &
19 ,dz8w,rho_phy &
20 ,alb,t3d,qv3d,qc3d,qr3d &
21 ,qi3d,qs3d,qg3d,qndrop3d &
22 ,p3d,p8w3d,pi3d,cldfra3d,rswtoa &
23 ,gmt,cp,g,julday,xtime,declin,solcon &
24 ,radfrq,degrad,taucldi,taucldc,warm_rain &
25 ,tauaer300,tauaer400,tauaer600,tauaer999 & ! jcb
26 ,gaer300,gaer400,gaer600,gaer999 & ! jcb
27 ,waer300,waer400,waer600,waer999 & ! jcb
28 ,aer_ra_feedback &
29 ,f_qv,f_qc,f_qr,f_qi,f_qs,f_qg,f_qndrop &
30 ,ids,ide, jds,jde, kds,kde &
31 ,ims,ime, jms,jme, kms,kme &
32 ,its,ite, jts,jte, kts,kte &
33 ,cosz_urb2d,omg_urb2d ) !Optional urban
34 !------------------------------------------------------------------
35 IMPLICIT NONE
36 !------------------------------------------------------------------
37 INTEGER, PARAMETER :: np = 75
38
39 INTEGER, INTENT(IN ) :: ids,ide, jds,jde, kds,kde, &
40 ims,ime, jms,jme, kms,kme, &
41 its,ite, jts,jte, kts,kte
42 LOGICAL, INTENT(IN ) :: warm_rain
43
44 INTEGER, INTENT(IN ) :: JULDAY
45
46
47 REAL, INTENT(IN ) :: RADFRQ,DEGRAD, &
48 XTIME,DECLIN,SOLCON
49 !
50 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
51 INTENT(IN ) :: P3D, &
52 P8W3D, &
53 pi3D, &
54 T3D, &
55 dz8w, &
56 rho_phy, &
57 CLDFRA3D
58
59
60 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
61 INTENT(INOUT) :: RTHRATEN
62 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
63 OPTIONAL, &
64 INTENT(INOUT) :: taucldi, &
65 taucldc
66 !
67 REAL, DIMENSION( ims:ime, jms:jme ), &
68 INTENT(IN ) :: XLAT, &
69 XLONG, &
70 ALB
71 !
72 REAL, DIMENSION( ims:ime, jms:jme ), &
73 INTENT(INOUT) :: GSW, &
74 RSWTOA
75 !
76 REAL, INTENT(IN ) :: GMT,CP,G
77 !
78
79 !
80 ! Optional
81 !
82 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), OPTIONAL , &
83 INTENT(IN ) :: tauaer300,tauaer400,tauaer600,tauaer999, & ! jcb
84 gaer300,gaer400,gaer600,gaer999, & ! jcb
85 waer300,waer400,waer600,waer999 ! jcb
86
87 INTEGER, INTENT(IN ), OPTIONAL :: aer_ra_feedback
88
89 REAL, DIMENSION( ims:ime, kms:kme, jms:jme ), &
90 OPTIONAL, &
91 INTENT(IN ) :: &
92 QV3D, &
93 QC3D, &
94 QR3D, &
95 QI3D, &
96 QS3D, &
97 QG3D, &
98 QNDROP3D
99
100 LOGICAL, OPTIONAL, INTENT(IN ) :: &
101 F_QV,F_QC,F_QR,F_QI,F_QS,F_QG, &
102 F_QNDROP
103
104 ! LOCAL VARS
105
106 REAL, DIMENSION( its:ite ) :: &
107 ts, &
108 cosz, &
109 rsuvbm, &
110 rsuvdf, &
111 rsirbm, &
112 rsirdf, &
113 p400, &
114 p700
115
116 INTEGER, DIMENSION( its:ite ) :: &
117 ict, &
118 icb
119
120 REAL, DIMENSION( its:ite, kts-1:kte, 2 ) :: taucld
121
122 REAL, DIMENSION( its:ite, kts-1:kte+1 ) :: flx, &
123 flxd
124 !
125 REAL, DIMENSION( its:ite, kts-1:kte ) :: O3
126 !
127 REAL, DIMENSION( its:ite, kts-1:kte, 11 ) :: &
128 taual, &
129 ssaal, &
130 asyal
131
132 REAL, DIMENSION( its:ite, kts-1:kte, 2 ) :: &
133 reff, &
134 cwc
135 REAL, DIMENSION( its: ite, kts-1:kte+1 ) :: &
136 P8W2D
137 REAL, DIMENSION( its: ite, kts-1:kte ) :: &
138 TTEN2D, &
139 qndrop2d, &
140 SH2D, &
141 P2D, &
142 T2D, &
143 fcld2D
144 real, DIMENSION( its:ite , kts:kte+1 ) :: phyd
145 real, DIMENSION( its:ite , kts:kte ) :: phydmid
146
147 REAL, DIMENSION( np, 5 ) :: pres, &
148 ozone
149 REAL, DIMENSION( np ) :: p
150
151 LOGICAL :: cldwater,overcast, predicate
152 !
153 INTEGER :: i,j,K,NK,ib,kk,mix,mkx
154
155 ! iprof = 1 : mid-latitude summer profile
156 ! = 2 : mid-latitude winter profile
157 ! = 3 : sub-arctic summer profile
158 ! = 4 : sub-arctic winter profile
159 ! = 5 : tropical profile
160 !
161
162 INTEGER :: iprof, &
163 is_summer, &
164 ie_summer, &
165 lattmp
166
167
168 !
169 REAL :: XLAT0,XLONG0
170 REAL :: fac,latrmp
171 REAL :: xt24,tloctm,hrang,xxlat
172
173 !URBAN
174 REAL, OPTIONAL, DIMENSION( ims:ime, jms:jme ), INTENT(OUT) :: COSZ_URB2D !urban
175 REAL, OPTIONAL, DIMENSION( ims:ime, jms:jme ), INTENT(OUT) :: OMG_URB2D !urban
176
177 real, dimension(11) :: midbands ! jcb
178 data midbands/.2,.235,.27,.2875,.3025,.305,.3625,.55,1.92,1.745,6.135/ ! jcb
179 real :: ang,slope ! jcb
180 character(len=200) :: msg !wig
181 real pi, third, relconst, lwpmin, rhoh2o
182 !
183 !--------------------------------------------------------------------------------
184 ! data set 1
185 ! mid-latitude summer (75 levels) : p(mb) o3(g/g)
186 ! surface temp = 294.0
187 !
188 data (pres(i,1),i=1,np)/ &
189 0.0006244, 0.0008759, 0.0012286, 0.0017234, 0.0024174, &
190 0.0033909, 0.0047565, 0.0066720, 0.0093589, 0.0131278, &
191 0.0184145, 0.0258302, 0.0362323, 0.0508234, 0.0712906, &
192 0.1000000, 0.1402710, 0.1967600, 0.2759970, 0.3871430, &
193 0.5430, 0.7617, 1.0685, 1.4988, 2.1024, 2.9490, &
194 4.1366, 5.8025, 8.1392, 11.4170, 16.0147, 22.4640, &
195 31.5105, 44.2001, 62.0000, 85.7750, 109.5500, 133.3250, &
196 157.1000, 180.8750, 204.6500, 228.4250, 252.2000, 275.9750, &
197 299.7500, 323.5250, 347.3000, 371.0750, 394.8500, 418.6250, &
198 442.4000, 466.1750, 489.9500, 513.7250, 537.5000, 561.2750, &
199 585.0500, 608.8250, 632.6000, 656.3750, 680.1500, 703.9250, &
200 727.7000, 751.4750, 775.2500, 799.0250, 822.8000, 846.5750, &
201 870.3500, 894.1250, 917.9000, 941.6750, 965.4500, 989.2250, &
202 1013.0000/
203 !
204 data (ozone(i,1),i=1,np)/ &
205 0.1793E-06, 0.2228E-06, 0.2665E-06, 0.3104E-06, 0.3545E-06, &
206 0.3989E-06, 0.4435E-06, 0.4883E-06, 0.5333E-06, 0.5786E-06, &
207 0.6241E-06, 0.6698E-06, 0.7157E-06, 0.7622E-06, 0.8557E-06, &
208 0.1150E-05, 0.1462E-05, 0.1793E-05, 0.2143E-05, 0.2512E-05, &
209 0.2902E-05, 0.3313E-05, 0.4016E-05, 0.5193E-05, 0.6698E-05, &
210 0.8483E-05, 0.9378E-05, 0.9792E-05, 0.1002E-04, 0.1014E-04, &
211 0.9312E-05, 0.7834E-05, 0.6448E-05, 0.5159E-05, 0.3390E-05, &
212 0.1937E-05, 0.1205E-05, 0.8778E-06, 0.6935E-06, 0.5112E-06, &
213 0.3877E-06, 0.3262E-06, 0.2770E-06, 0.2266E-06, 0.2020E-06, &
214 0.1845E-06, 0.1679E-06, 0.1519E-06, 0.1415E-06, 0.1317E-06, &
215 0.1225E-06, 0.1137E-06, 0.1055E-06, 0.1001E-06, 0.9487E-07, &
216 0.9016E-07, 0.8641E-07, 0.8276E-07, 0.7930E-07, 0.7635E-07, &
217 0.7347E-07, 0.7065E-07, 0.6821E-07, 0.6593E-07, 0.6368E-07, &
218 0.6148E-07, 0.5998E-07, 0.5859E-07, 0.5720E-07, 0.5582E-07, &
219 0.5457E-07, 0.5339E-07, 0.5224E-07, 0.5110E-07, 0.4999E-07/
220
221 !--------------------------------------------------------------------------------
222 ! data set 2
223 ! mid-latitude winter (75 levels) : p(mb) o3(g/g)
224 ! surface temp = 272.2
225 !
226 data (pres(i,2),i=1,np)/ &
227 0.0006244, 0.0008759, 0.0012286, 0.0017234, 0.0024174, &
228 0.0033909, 0.0047565, 0.0066720, 0.0093589, 0.0131278, &
229 0.0184145, 0.0258302, 0.0362323, 0.0508234, 0.0712906, &
230 0.1000000, 0.1402710, 0.1967600, 0.2759970, 0.3871430, &
231 0.5430, 0.7617, 1.0685, 1.4988, 2.1024, 2.9490, &
232 4.1366, 5.8025, 8.1392, 11.4170, 16.0147, 22.4640, &
233 31.5105, 44.2001, 62.0000, 85.9000, 109.8000, 133.7000, &
234 157.6000, 181.5000, 205.4000, 229.3000, 253.2000, 277.1000, &
235 301.0000, 324.9000, 348.8000, 372.7000, 396.6000, 420.5000, &
236 444.4000, 468.3000, 492.2000, 516.1000, 540.0000, 563.9000, &
237 587.8000, 611.7000, 635.6000, 659.5000, 683.4000, 707.3000, &
238 731.2000, 755.1000, 779.0000, 802.9000, 826.8000, 850.7000, &
239 874.6000, 898.5000, 922.4000, 946.3000, 970.2000, 994.1000, &
240 1018.0000/
241 !
242 data (ozone(i,2),i=1,np)/ &
243 0.2353E-06, 0.3054E-06, 0.3771E-06, 0.4498E-06, 0.5236E-06, &
244 0.5984E-06, 0.6742E-06, 0.7511E-06, 0.8290E-06, 0.9080E-06, &
245 0.9881E-06, 0.1069E-05, 0.1152E-05, 0.1319E-05, 0.1725E-05, &
246 0.2145E-05, 0.2581E-05, 0.3031E-05, 0.3497E-05, 0.3980E-05, &
247 0.4478E-05, 0.5300E-05, 0.6725E-05, 0.8415E-05, 0.1035E-04, &
248 0.1141E-04, 0.1155E-04, 0.1143E-04, 0.1093E-04, 0.1060E-04, &
249 0.9720E-05, 0.8849E-05, 0.7424E-05, 0.6023E-05, 0.4310E-05, &
250 0.2820E-05, 0.1990E-05, 0.1518E-05, 0.1206E-05, 0.9370E-06, &
251 0.7177E-06, 0.5450E-06, 0.4131E-06, 0.3277E-06, 0.2563E-06, &
252 0.2120E-06, 0.1711E-06, 0.1524E-06, 0.1344E-06, 0.1199E-06, &
253 0.1066E-06, 0.9516E-07, 0.8858E-07, 0.8219E-07, 0.7598E-07, &
254 0.6992E-07, 0.6403E-07, 0.5887E-07, 0.5712E-07, 0.5540E-07, &
255 0.5370E-07, 0.5214E-07, 0.5069E-07, 0.4926E-07, 0.4785E-07, &
256 0.4713E-07, 0.4694E-07, 0.4676E-07, 0.4658E-07, 0.4641E-07, &
257 0.4634E-07, 0.4627E-07, 0.4619E-07, 0.4612E-07, 0.4605E-07/
258
259
260 !--------------------------------------------------------------------------------
261 ! data set 3
262 ! sub-arctic summer (75 levels) : p(mb) o3(g/g)
263 ! surface temp = 287.0
264 !
265 data (pres(i,3),i=1,np)/ &
266 0.0006244, 0.0008759, 0.0012286, 0.0017234, 0.0024174, &
267 0.0033909, 0.0047565, 0.0066720, 0.0093589, 0.0131278, &
268 0.0184145, 0.0258302, 0.0362323, 0.0508234, 0.0712906, &
269 0.1000000, 0.1402710, 0.1967600, 0.2759970, 0.3871430, &
270 0.5430, 0.7617, 1.0685, 1.4988, 2.1024, 2.9490, &
271 4.1366, 5.8025, 8.1392, 11.4170, 16.0147, 22.4640, &
272 31.5105, 44.2001, 62.0000, 85.7000, 109.4000, 133.1000, &
273 156.8000, 180.5000, 204.2000, 227.9000, 251.6000, 275.3000, &
274 299.0000, 322.7000, 346.4000, 370.1000, 393.8000, 417.5000, &
275 441.2000, 464.9000, 488.6000, 512.3000, 536.0000, 559.7000, &
276 583.4000, 607.1000, 630.8000, 654.5000, 678.2000, 701.9000, &
277 725.6000, 749.3000, 773.0000, 796.7000, 820.4000, 844.1000, &
278 867.8000, 891.5000, 915.2000, 938.9000, 962.6000, 986.3000, &
279 1010.0000/
280 !
281 data (ozone(i,3),i=1,np)/ &
282 0.1728E-06, 0.2131E-06, 0.2537E-06, 0.2944E-06, 0.3353E-06, &
283 0.3764E-06, 0.4176E-06, 0.4590E-06, 0.5006E-06, 0.5423E-06, &
284 0.5842E-06, 0.6263E-06, 0.6685E-06, 0.7112E-06, 0.7631E-06, &
285 0.1040E-05, 0.1340E-05, 0.1660E-05, 0.2001E-05, 0.2362E-05, &
286 0.2746E-05, 0.3153E-05, 0.3762E-05, 0.4988E-05, 0.6518E-05, &
287 0.8352E-05, 0.9328E-05, 0.9731E-05, 0.8985E-05, 0.7632E-05, &
288 0.6814E-05, 0.6384E-05, 0.5718E-05, 0.4728E-05, 0.4136E-05, &
289 0.3033E-05, 0.2000E-05, 0.1486E-05, 0.1121E-05, 0.8680E-06, &
290 0.6474E-06, 0.5164E-06, 0.3921E-06, 0.2996E-06, 0.2562E-06, &
291 0.2139E-06, 0.1723E-06, 0.1460E-06, 0.1360E-06, 0.1267E-06, &
292 0.1189E-06, 0.1114E-06, 0.1040E-06, 0.9678E-07, 0.8969E-07, &
293 0.8468E-07, 0.8025E-07, 0.7590E-07, 0.7250E-07, 0.6969E-07, &
294 0.6694E-07, 0.6429E-07, 0.6208E-07, 0.5991E-07, 0.5778E-07, &
295 0.5575E-07, 0.5403E-07, 0.5233E-07, 0.5067E-07, 0.4904E-07, &
296 0.4721E-07, 0.4535E-07, 0.4353E-07, 0.4173E-07, 0.3997E-07/
297
298
299 !--------------------------------------------------------------------------------
300 ! data set 3
301 ! sub-arctic winter (75 levels) : p(mb) o3(g/g)
302 ! surface temp = 257.1
303 !
304 data (pres(i,4),i=1,np)/ &
305 0.0006244, 0.0008759, 0.0012286, 0.0017234, 0.0024174, &
306 0.0033909, 0.0047565, 0.0066720, 0.0093589, 0.0131278, &
307 0.0184145, 0.0258302, 0.0362323, 0.0508234, 0.0712906, &
308 0.1000000, 0.1402710, 0.1967600, 0.2759970, 0.3871430, &
309 0.5430, 0.7617, 1.0685, 1.4988, 2.1024, 2.9490, &
310 4.1366, 5.8025, 8.1392, 11.4170, 16.0147, 22.4640, &
311 31.5105, 44.2001, 62.0000, 85.7750, 109.5500, 133.3250, &
312 157.1000, 180.8750, 204.6500, 228.4250, 252.2000, 275.9750, &
313 299.7500, 323.5250, 347.3000, 371.0750, 394.8500, 418.6250, &
314 442.4000, 466.1750, 489.9500, 513.7250, 537.5000, 561.2750, &
315 585.0500, 608.8250, 632.6000, 656.3750, 680.1500, 703.9250, &
316 727.7000, 751.4750, 775.2500, 799.0250, 822.8000, 846.5750, &
317 870.3500, 894.1250, 917.9000, 941.6750, 965.4500, 989.2250, &
318 1013.0000/
319 !
320 data (ozone(i,4),i=1,np)/ &
321 0.2683E-06, 0.3562E-06, 0.4464E-06, 0.5387E-06, 0.6333E-06, &
322 0.7301E-06, 0.8291E-06, 0.9306E-06, 0.1034E-05, 0.1140E-05, &
323 0.1249E-05, 0.1360E-05, 0.1474E-05, 0.1855E-05, 0.2357E-05, &
324 0.2866E-05, 0.3383E-05, 0.3906E-05, 0.4437E-05, 0.4975E-05, &
325 0.5513E-05, 0.6815E-05, 0.8157E-05, 0.1008E-04, 0.1200E-04, &
326 0.1242E-04, 0.1250E-04, 0.1157E-04, 0.1010E-04, 0.9063E-05, &
327 0.8836E-05, 0.8632E-05, 0.8391E-05, 0.7224E-05, 0.6054E-05, &
328 0.4503E-05, 0.3204E-05, 0.2278E-05, 0.1833E-05, 0.1433E-05, &
329 0.9996E-06, 0.7440E-06, 0.5471E-06, 0.3944E-06, 0.2852E-06, &
330 0.1977E-06, 0.1559E-06, 0.1333E-06, 0.1126E-06, 0.9441E-07, &
331 0.7678E-07, 0.7054E-07, 0.6684E-07, 0.6323E-07, 0.6028E-07, &
332 0.5746E-07, 0.5468E-07, 0.5227E-07, 0.5006E-07, 0.4789E-07, &
333 0.4576E-07, 0.4402E-07, 0.4230E-07, 0.4062E-07, 0.3897E-07, &
334 0.3793E-07, 0.3697E-07, 0.3602E-07, 0.3506E-07, 0.3413E-07, &
335 0.3326E-07, 0.3239E-07, 0.3153E-07, 0.3069E-07, 0.2987E-07/
336
337 !--------------------------------------------------------------------------------
338 ! data set 4
339 ! tropical (75 levels) : p(mb) o3(g/g)
340 ! surface temp = 300.0
341 !
342 data (pres(i,5),i=1,np)/ &
343 0.0006244, 0.0008759, 0.0012286, 0.0017234, 0.0024174, &
344 0.0033909, 0.0047565, 0.0066720, 0.0093589, 0.0131278, &
345 0.0184145, 0.0258302, 0.0362323, 0.0508234, 0.0712906, &
346 0.1000000, 0.1402710, 0.1967600, 0.2759970, 0.3871430, &
347 0.5430, 0.7617, 1.0685, 1.4988, 2.1024, 2.9490, &
348 4.1366, 5.8025, 8.1392, 11.4170, 16.0147, 22.4640, &
349 31.5105, 44.2001, 62.0000, 85.7750, 109.5500, 133.3250, &
350 157.1000, 180.8750, 204.6500, 228.4250, 252.2000, 275.9750, &
351 299.7500, 323.5250, 347.3000, 371.0750, 394.8500, 418.6250, &
352 442.4000, 466.1750, 489.9500, 513.7250, 537.5000, 561.2750, &
353 585.0500, 608.8250, 632.6000, 656.3750, 680.1500, 703.9250, &
354 727.7000, 751.4750, 775.2500, 799.0250, 822.8000, 846.5750, &
355 870.3500, 894.1250, 917.9000, 941.6750, 965.4500, 989.2250, &
356 1013.0000/
357 !
358 data (ozone(i,5),i=1,np)/ &
359 0.1993E-06, 0.2521E-06, 0.3051E-06, 0.3585E-06, 0.4121E-06, &
360 0.4661E-06, 0.5203E-06, 0.5748E-06, 0.6296E-06, 0.6847E-06, &
361 0.7402E-06, 0.7959E-06, 0.8519E-06, 0.9096E-06, 0.1125E-05, &
362 0.1450E-05, 0.1794E-05, 0.2156E-05, 0.2538E-05, 0.2939E-05, &
363 0.3362E-05, 0.3785E-05, 0.4753E-05, 0.6005E-05, 0.7804E-05, &
364 0.9635E-05, 0.1023E-04, 0.1067E-04, 0.1177E-04, 0.1290E-04, &
365 0.1134E-04, 0.9223E-05, 0.6667E-05, 0.3644E-05, 0.1545E-05, &
366 0.5355E-06, 0.2523E-06, 0.2062E-06, 0.1734E-06, 0.1548E-06, &
367 0.1360E-06, 0.1204E-06, 0.1074E-06, 0.9707E-07, 0.8960E-07, &
368 0.8419E-07, 0.7962E-07, 0.7542E-07, 0.7290E-07, 0.7109E-07, &
369 0.6940E-07, 0.6786E-07, 0.6635E-07, 0.6500E-07, 0.6370E-07, &
370 0.6244E-07, 0.6132E-07, 0.6022E-07, 0.5914E-07, 0.5884E-07, &
371 0.5855E-07, 0.5823E-07, 0.5772E-07, 0.5703E-07, 0.5635E-07, &
372 0.5570E-07, 0.5492E-07, 0.5412E-07, 0.5335E-07, 0.5260E-07, &
373 0.5167E-07, 0.5063E-07, 0.4961E-07, 0.4860E-07, 0.4761E-07/
374
375 !--------------------------------------------------------------------------------
376
377 #ifdef WRF_CHEM
378 IF ( aer_ra_feedback == 1) then
379 IF ( .NOT. &
380 ( PRESENT(tauaer300) .AND. &
381 PRESENT(tauaer400) .AND. &
382 PRESENT(tauaer600) .AND. &
383 PRESENT(tauaer999) .AND. &
384 PRESENT(gaer300) .AND. &
385 PRESENT(gaer400) .AND. &
386 PRESENT(gaer600) .AND. &
387 PRESENT(gaer999) .AND. &
388 PRESENT(waer300) .AND. &
389 PRESENT(waer400) .AND. &
390 PRESENT(waer600) .AND. &
391 PRESENT(waer999) ) ) THEN
392 CALL wrf_error_fatal ( 'Warning: missing fields required for aerosol radiation' )
393 ENDIF
394 ENDIF
395 #endif
396 cldwater = .true.
397 overcast = .false.
398
399 mix=ite-its+1
400 mkx=kte-kts+1
401
402 is_summer=80
403 ie_summer=265
404
405 ! testing, need to change iprof, which is function of lat and julian day
406 ! iprof = 1 : mid-latitude summer profile
407 ! = 2 : mid-latitude winter profile
408 ! = 3 : sub-arctic summer profile
409 ! = 4 : sub-arctic winter profile
410 ! = 5 : tropical profile
411
412 IF (abs(center_lat) .le. 30. ) THEN ! tropic
413 iprof = 5
414 ELSE
415 IF (center_lat .gt. 0.) THEN
416 IF (center_lat .gt. 60. ) THEN ! arctic
417 IF (JULDAY .gt. is_summer .and. JULDAY .lt. ie_summer ) THEN
418 ! arctic summer
419 iprof = 3
420 ELSE
421 ! arctic winter
422 iprof = 4
423 ENDIF
424 ELSE ! midlatitude
425 IF (JULDAY .gt. is_summer .and. JULDAY .lt. ie_summer ) THEN
426 ! north midlatitude summer
427 iprof = 1
428 ELSE
429 ! north midlatitude winter
430 iprof = 2
431 ENDIF
432 ENDIF
433
434 ELSE
435 IF (center_lat .lt. -60. ) THEN ! antarctic
436 IF (JULDAY .lt. is_summer .or. JULDAY .gt. ie_summer ) THEN
437 ! antarctic summer
438 iprof = 3
439 ELSE
440 ! antarctic winter
441 iprof = 4
442 ENDIF
443 ELSE ! midlatitude
444 IF (JULDAY .lt. is_summer .or. JULDAY .gt. ie_summer ) THEN
445 ! south midlatitude summer
446 iprof = 1
447 ELSE
448 ! south midlatitude winter
449 iprof = 2
450 ENDIF
451 ENDIF
452
453 ENDIF
454 ENDIF
455
456
457 j_loop: DO J=jts,jte
458
459 DO K=kts,kte
460 DO I=its,ite
461 cwc(i,k,1) = 0.
462 cwc(i,k,2) = 0.
463 ENDDO
464 ENDDO
465
466 DO K=1,np
467 p(k)=pres(k,iprof)
468 ENDDO
469
470 do k = kts,kte+1
471 do i = its,ite
472 if(k.eq.kts)then
473 phyd(i,k)=p8w3d(i,kts,j)
474 else
475 phyd(i,k)=phyd(i,k-1) - g*rho_phy(i,k-1,j)*dz8w(i,k-1,j)
476 phydmid(i,k-1)=0.5*(phyd(i,k-1)+phyd(i,k))
477 endif
478 enddo
479 enddo
480
481 ! reverse vars
482 !
483 DO K=kts,kte+1
484 DO I=its,ite
485 NK=kme-K+kms
486 P8W2D(I,K)=phyd(I,NK)*0.01 ! P8w2D is in mb
487 ENDDO
488 ENDDO
489
490 DO I=its,ite
491 P8W2D(I,0)=.0
492 ENDDO
493 !
494 DO K=kts,kte
495 DO I=its,ite
496 NK=kme-1-K+kms
497 TTEN2D(I,K)=0.
498 T2D(I,K)=T3D(I,NK,J)
499
500 ! SH2D specific humidity
501 SH2D(I,K)=QV3D(I,NK,J)/(1.+QV3D(I,NK,J))
502 SH2D(I,K)=max(0.,SH2D(I,K))
503 cwc(I,K,2)=QC3D(I,NK,J)
504 cwc(I,K,2)=max(0.,cwc(I,K,2))
505
506 P2D(I,K)=phydmid(I,NK)*0.01 ! P2D is in mb
507 fcld2D(I,K)=CLDFRA3D(I,NK,J)
508 ENDDO
509 ENDDO
510
511 ! This logic is tortured because cannot test F_QI unless
512 ! it is present, and order of evaluation of expressions
513 ! is not specified in Fortran
514
515 IF ( PRESENT ( F_QI ) ) THEN
516 predicate = F_QI
517 ELSE
518 predicate = .FALSE.
519 ENDIF
520
521 IF (.NOT. warm_rain .AND. .NOT. predicate ) THEN
522 DO K=kts,kte
523 DO I=its,ite
524 IF (T2D(I,K) .lt. 273.15) THEN
525 cwc(I,K,1)=cwc(I,K,2)
526 cwc(I,K,2)=0.
527 ENDIF
528 ENDDO
529 ENDDO
530 ENDIF
531
532 IF ( PRESENT( F_QNDROP ) ) THEN
533 IF ( F_QNDROP ) THEN
534 DO K=kts,kte
535 DO I=its,ite
536 NK=kme-1-K+kms
537 qndrop2d(I,K)=qndrop3d(I,NK,j)
538 ENDDO
539 ENDDO
540 qndrop2d(:,kts-1)=0.
541 END IF
542 END IF
543
544 DO I=its,ite
545 TTEN2D(I,0)=0.
546 T2D(I,0)=T2D(I,1)
547 ! SH2D specific humidity
548 SH2D(I,0)=0.5*SH2D(i,1)
549 cwc(I,0,2)=0.
550 cwc(I,0,1)=0.
551 P2D(I,0)=0.5*(P8W2D(I,0)+P8W2D(I,1))
552 fcld2D(I,0)=0.
553 ENDDO
554 !
555 IF ( PRESENT( F_QI ) .AND. PRESENT( qi3d) ) THEN
556 IF ( (F_QI) ) THEN
557 DO K=kts,kte
558 DO I=its,ite
559 NK=kme-1-K+kms
560 cwc(I,K,1)=QI3D(I,NK,J)
561 cwc(I,K,1)=max(0.,cwc(I,K,1))
562 ENDDO
563 ENDDO
564 ENDIF
565 ENDIF
566 !
567 ! ... Vertical profiles for ozone
568 !
569 call o3prof (np, p, ozone(1,iprof), its, ite, kts-1, kte, P2D, O3)
570
571 ! ... Vertical profiles for effective particle size
572 !
573 pi = 4.*atan(1.0)
574 third=1./3.
575 rhoh2o=1.e3
576 relconst=3/(4.*pi*rhoh2o)
577 ! minimun liquid water path to calculate rel
578 ! corresponds to optical depth of 1.e-3 for radius 4 microns.
579 lwpmin=3.e-5
580 do k = kts-1, kte
581 do i = its, ite
582 reff(i,k,2) = 10.
583 if( PRESENT( F_QNDROP ) ) then
584 if( F_QNDROP ) then
585 if ( cwc(i,k,2)*(P8W2D(I,K+1)-P8W2D(I,K)).gt.lwpmin.and. &
586 qndrop2d(i,k).gt.1000. ) then
587 reff(i,k,2)=(relconst*cwc(i,k,2)/qndrop2d(i,k))**third ! effective radius in m
588 ! apply scaling from Martin et al., JAS 51, 1830.
589 reff(i,k,2)=1.1*reff(i,k,2)
590 reff(i,k,2)=reff(i,k,2)*1.e6 ! convert from m to microns
591 reff(i,k,2)=max(reff(i,k,2),4.)
592 reff(i,k,2)=min(reff(i,k,2),20.)
593 end if
594 end if
595 end if
596 reff(i,k,1) = 80.
597 end do
598 end do
599 !
600 ! ... Level indices separating high, middle and low clouds
601 !
602 do i = its, ite
603 p400(i) = 1.e5
604 p700(i) = 1.e5
605 enddo
606
607 do k = kts-1,kte+1
608 do i = its, ite
609 if (abs(P8W2D(i,k) - 400.) .lt. p400(i)) then
610 p400(i) = abs(P8W2D(i,k) - 400.)
611 ict(i) = k
612 endif
613 if (abs(P8W2D(i,k) - 700.) .lt. p700(i)) then
614 p700(i) = abs(P8W2D(i,k) - 700.)
615 icb(i) = k
616 endif
617 end do
618 end do
619
620 !wig beg
621 ! ... Aerosol effects. Added aerosol feedbacks with MOSAIC, Dec. 2005.
622 !
623 do ib = 1, 11
624 do k = kts-1,kte
625 do i = its,ite
626 taual(i,k,ib) = 0.
627 ssaal(i,k,ib) = 0.
628 asyal(i,k,ib) = 0.
629 end do
630 end do
631 end do
632
633 #ifdef WRF_CHEM
634 IF ( AER_RA_FEEDBACK == 1) then
635 !wig end
636 do ib = 1, 11
637 do k = kts-1,kte-1 !wig
638 do i = its,ite
639
640 ! taual(i,kte-k,ib) = 0.
641 ! ssaal(i,kte-k,ib) = 0.
642 ! asyal(i,kte-k,ib) = 0.
643
644 !jcb beg
645 ! convert optical properties at 300,400,600, and 999 to conform to the band wavelengths
646 ! these are: 200,235,270,287.5,302.5,305,362.5,550,1920,1745,6135; why the emphasis on the UV?
647 ! taual - use angstrom exponent
648 if(tauaer300(i,k+1,j).gt.thresh .and. tauaer999(i,k+1,j).gt.thresh) then
649 ang=log(tauaer300(i,k+1,j)/tauaer999(i,k+1,j))/log(999./300.)
650 ! write(6,*)i,k,ang,tauaer300(i,k+1,j),tauaer999(i,k+1,j)
651 taual(i,kte-k,ib)=tauaer400(i,k+1,j)*(0.4/midbands(ib))**ang ! notice reserved variable
652 ! write(6,10001)i,k,ang,tauaer300(i,k+1,j),tauaer999(i,k+1,j),midbands(ib),taual(i,k,ib)
653 !10001 format(i3,i3,5f12.6)
654
655 ! ssa - linear interpolation; extrapolation
656 slope=(waer600(i,k+1,j)-waer400(i,k+1,j))/.2
657 ssaal(i,kte-k,ib) = slope*(midbands(ib)-.6)+waer600(i,k+1,j) ! notice reversed variables
658 if(ssaal(i,kte-k,ib).lt.0.4) ssaal(i,kte-k,ib)=0.4
659 if(ssaal(i,kte-k,ib).ge.1.0) ssaal(i,kte-k,ib)=1.0
660
661 ! g - linear interpolation;extrapolation
662 slope=(gaer600(i,k+1,j)-gaer400(i,k+1,j))/.2
663 asyal(i,kte-k,ib) = slope*(midbands(ib)-.6)+gaer600(i,k+1,j) ! notice reversed varaibles
664 if(asyal(i,kte-k,ib).lt.0.5) asyal(i,kte-k,ib)=0.5
665 if(asyal(i,kte-k,ib).ge.1.0) asyal(i,kte-k,ib)=1.0
666 endif
667 !jcb end
668 end do
669 end do
670 end do
671
672 !wig beg
673 do ib = 1, 11
674 do i = its,ite
675 slope = 0. !use slope as a sum holder
676 do k = kts-1,kte
677 slope = slope + taual(i,k,ib)
678 end do
679 if( slope < 0. ) then
680 write(msg,'("ERROR: Negative total optical depth of ",f8.2," at point i,j,ib=",3i5)') slope,i,j,ib
681 call wrf_error_fatal(msg)
682 else if( slope > 5. ) then
683 call wrf_message("-------------------------")
684 write(msg,'("WARNING: Large total optical depth of ",f8.2," at point i,j,ib=",3i5)') slope,i,j,ib
685 call wrf_message(msg)
686
687 call wrf_message("Diagnostics 1: k, tauaer300, tauaer400, tauaer600, tauaer999")
688 do k=kts,kte
689 write(msg,'(i4,4f8.2)') k, tauaer300(i,k,j), tauaer400(i,k,j), &
690 tauaer600(i,k,j), tauaer999(i,k,j)
691 call wrf_message(msg)
692 end do
693
694 call wrf_message("Diagnostics 2: k, gaer300, gaer400, gaer600, gaer999")
695 do k=kts,kte
696 write(msg,'(i4,4f8.2)') k, gaer300(i,k,j), gaer400(i,k,j), &
697 gaer600(i,k,j), gaer999(i,k,j)
698 call wrf_message(msg)
699 end do
700
701 call wrf_message("Diagnostics 3: k, waer300, waer400, waer600, waer999")
702 do k=kts,kte
703 write(msg,'(i4,4f8.2)') k, waer300(i,k,j), waer400(i,k,j), &
704 waer600(i,k,j), waer999(i,k,j)
705 call wrf_message(msg)
706 end do
707
708 call wrf_message("Diagnostics 4: k, ssaal, asyal, taual")
709 do k=kts-1,kte
710 write(msg,'(i4,3f8.2)') k, ssaal(i,k,ib), asyal(i,k,ib), taual(i,k,ib)
711 call wrf_message(msg)
712 end do
713 call wrf_message("-------------------------")
714 end if
715 end do
716 end do
717 !wig end
718 endif
719 #endif
720 !
721 ! ... Initialize output arrays
722 !
723 do ib = 1, 2
724 do k = kts-1, kte
725 do i = its, ite
726 taucld(i,k,ib) = 0.
727 end do
728 end do
729 end do
730 !
731 do k = kts-1,kte+1
732 do i = its,ite
733 flx(i,k) = 0.
734 flxd(i,k) = 0.
735 end do
736 end do
737 !
738 ! ... Solar zenith angle
739 !
740 do i = its,ite
741 xt24 = mod(xtime + radfrq * 0.5, 1440.)
742 tloctm = GMT + xt24 / 60. + XLONG(i,j) / 15.
743 hrang = 15. * (tloctm - 12.) * degrad
744 xxlat = XLAT(i,j) * degrad
745 cosz(i) = sin(xxlat) * sin(declin) + &
746 cos(xxlat) * cos(declin) * cos(hrang)
747 !urban
748 if(present(COSZ_URB2D)) COSZ_URB2D(i,j)=cosz(i) !urban
749 if(present(OMG_URB2D)) OMG_URB2D(i,j)=hrang !urban
750 rsuvbm(i) = ALB(i,j)
751 rsuvdf(i) = ALB(i,j)
752 rsirbm(i) = ALB(i,j)
753 rsirdf(i) = ALB(i,j)
754 end do
755
756 call sorad (mix,1,1,mkx+1,p8w2D,t2D,sh2D,o3, &
757 overcast,cldwater,cwc,taucld,reff,fcld2D,ict,icb,&
758 taual,ssaal,asyal, &
759 cosz,rsuvbm,rsuvdf,rsirbm,rsirdf, &
760 flx,flxd)
761 !
762 ! ... Convert the units of flx and flc from fraction to w/m^2
763 !
764 do k = kts, kte
765 do i = its, ite
766 nk=kme-1-k+kms
767 if(present(taucldc)) taucldc(i,nk,j)=taucld(i,k,2)
768 if(present(taucldi)) taucldi(i,nk,j)=taucld(i,k,1)
769 enddo
770 enddo
771
772 do k = kts, kte+1
773 do i = its, ite
774 if (cosz(i) .lt. thresh) then
775 flx(i,k) = 0.
776 else
777 flx(i,k) = flx(i,k) * SOLCON * cosz(i)
778 endif
779 end do
780 end do
781 !
782 ! ... Calculate heating rate (deg/sec)
783 !
784 fac = .01 * g / Cp
785 do k = kts, kte
786 do i = its, ite
787 if (cosz(i) .gt. thresh) then
788 TTEN2D(i,k) = - fac * (flx(i,k) - flx(i,k+1))/ &
789 (p8w2d(i,k)-p8w2d(i,k+1))
790 endif
791 end do
792 end do
793
794 ! upward top of atmosphere
795 do i = its, ite
796 if (cosz(i) .le. thresh) then
797 RSWTOA(i,j) = 0.
798 else
799 RSWTOA(i,j) = flx(i,kts) - flxd(i,kts) * SOLCON * cosz(i)
800 ! print *,'cosz,rswtoa=',cosz(i),rswtoa(i,j)
801 endif
802 end do
803 !
804 ! ... Absorbed part in surface energy budget
805 !
806 do i = its, ite
807 if (cosz(i) .le. thresh) then
808 GSW(i,j) = 0.
809 else
810 GSW(i,j) = (1. - rsuvbm(i)) * flxd(i,kte+1) * SOLCON * cosz(i)
811 endif
812 end do
813
814 DO K=kts,kte
815 NK=kme-1-K+kms
816 DO I=its,ite
817 RTHRATEN(I,K,J)=RTHRATEN(I,K,J)+TTEN2D(I,NK)/pi3D(I,K,J)
818 ENDDO
819 ENDDO
820 !
821 ENDDO j_loop
822
823 END SUBROUTINE GSFCSWRAD
824
825 !********************* Version Solar-6 (May 8, 1997) *****************
826
827 subroutine sorad (m,n,ndim,np,pl,ta,wa,oa, &
828 overcast,cldwater,cwc,taucld,reff,fcld,ict,icb, &
829 taual,ssaal,asyal, &
830 cosz,rsuvbm,rsuvdf,rsirbm,rsirdf, &
831 flx,flxd)
832
833 !************************************************************************
834 !
835 ! Version Solar-6 (May 8, 1997)
836 !
837 ! New feature of this version is:
838 ! (1) An option is added for scaling the cloud optical thickness. If
839 ! the fractional cloud cover, fcld, in an atmospheric model is alway
840 ! either 1 or 0 (i.e. partly cloudy sky is not allowed), it does
841 ! not require the scaling of cloud optical thickness, and the
842 ! option "overcast" can be set to .true. Computation is faster
843 ! with this option than with overcast=.false.
844 !
845 !**********************************************************************
846 !
847 ! Version Solar-5 (April 1997)
848 !
849 ! New features of this version are:
850 ! (1) Cloud optical properties can be computed from cloud water/ice
851 ! amount and the effective particle size.
852 ! (2) Aerosol optical properties are functions of height and band.
853 ! (3) A maximum-random cloud overlapping approximation is applied.
854 !
855 !*********************************************************************
856 !
857 ! This routine computes solar fluxes due to the absoption by water
858 ! vapor, ozone, co2, o2, clouds, and aerosols and due to the
859 ! scattering by clouds, aerosols, and gases.
860 !
861 ! The solar spectrum is divided into one UV+visible band and three IR
862 ! bands separated by the wavelength 0.7 micron. The UV+visible band
863 ! is further divided into eight sub-bands.
864 !
865 ! This is a vectorized code. It computes fluxes simultaneously for
866 ! (m x n) soundings, which is a subset of (m x ndim) soundings.
867 ! In a global climate model, m and ndim correspond to the numbers of
868 ! grid boxes in the zonal and meridional directions, respectively.
869 !
870 ! Ice and liquid cloud particles are allowed to co-exist in a layer.
871 !
872 ! There is an option of providing either cloud ice/water mixing ratio
873 ! (cwc) or thickness (taucld). If the former is provided, set
874 ! cldwater=.true., and taucld will be computed from cwc and reff as a
875 ! function of spectra band. Otherwise, set cldwater=.false., and
876 ! specify taucld, independent of spectral band.
877 !
878 ! If no information is available for reff, a default value of
879 ! 10 micron for liquid water and 75 micron for ice can be used.
880 ! For a clear layer, reff can be set to any values except zero.
881 !
882 ! The maximum-random assumption is applied for treating cloud
883 ! overlapping.
884
885 ! Clouds are grouped into high, middle, and low clouds separated by
886 ! the level indices ict and icb. For detail, see subroutine cldscale.
887 !
888 ! In a high spatial-resolution atmospheric model, fractional cloud cover
889 ! might be computed to be either 0 or 1. In such a case, scaling of the
890 ! cloud optical thickness is not necessary, and the computation can be
891 ! made faster by setting overcast=.true. The option overcast=.false.
892 ! can be applied to any values of the fractional cloud cover, but the
893 ! computation is slower.
894 !
895 ! Aerosol optical thickness, single-scattering albaedo, and asymmtry
896 ! factor can be specified as functions of height and spectral band.
897 !
898 !----- Input parameters:
899 ! units size
900 ! number of soundings in zonal direction (m) n/d 1
901 ! number of soundings in meridional direction (n) n/d 1
902 ! maximum number of soundings in n/d 1
903 ! meridional direction (ndim>=n)
904 ! number of atmospheric layers (np) n/d 1
905 ! level pressure (pl) mb m*ndim*(np+1)
906 ! layer temperature (ta) k m*ndim*np
907 ! layer specific humidity (wa) gm/gm m*ndim*np
908 ! layer ozone concentration (oa) gm/gm m*ndim*np
909 ! co2 mixing ratio by volumn (co2) pppv 1
910 ! option for scaling cloud optical thickness n/d 1
911 ! overcast="true" if scaling is NOT required
912 ! overcast="fasle" if scaling is required
913 ! option for cloud optical thickness n/d 1
914 ! cldwater="true" if cwc is provided
915 ! cldwater="false" if taucld is provided
916 ! cloud water mixing ratio (cwc) gm/gm m*ndim*np*2
917 ! index 1 for ice particles
918 ! index 2 for liquid drops
919 ! cloud optical thickness (taucld) n/d m*ndim*np*2
920 ! index 1 for ice particles
921 ! index 2 for liquid drops
922 ! effective cloud-particle size (reff) micrometer m*ndim*np*2
923 ! index 1 for ice particles
924 ! index 2 for liquid drops
925 ! cloud amount (fcld) fraction m*ndim*np
926 ! level index separating high and middle n/d 1
927 ! clouds (ict)
928 ! level index separating middle and low n/d 1
929 ! clouds (icb)
930 ! aerosol optical thickness (taual) n/d m*ndim*np*11
931 ! aerosol single-scattering albedo (ssaal) n/d m*ndim*np*11
932 ! aerosol asymmetry factor (asyal) n/d m*ndim*np*11
933 ! in the uv region :
934 ! index 1 for the 0.175-0.225 micron band
935 ! index 2 for the 0.225-0.245; 0.260-0.280 micron band
936 ! index 3 for the 0.245-0.260 micron band
937 ! index 4 for the 0.280-0.295 micron band
938 ! index 5 for the 0.295-0.310 micron band
939 ! index 6 for the 0.310-0.320 micron band
940 ! index 7 for the 0.325-0.400 micron band
941 ! in the par region :
942 ! index 8 for the 0.400-0.700 micron band
943 ! in the infrared region :
944 ! index 9 for the 0.700-1.220 micron band
945 ! index 10 for the 1.220-2.270 micron band
946 ! index 11 for the 2.270-10.00 micron band
947 ! cosine of solar zenith angle (cosz) n/d m*ndim
948 ! uv+visible sfc albedo for beam radiation
949 ! for wavelengths<0.7 micron (rsuvbm) fraction m*ndim
950 ! uv+visible sfc albedo for diffuse radiation
951 ! for wavelengths<0.7 micron (rsuvdf) fraction m*ndim
952 ! ir sfc albedo for beam radiation
953 ! for wavelengths>0.7 micron (rsirbm) fraction m*ndim
954 ! ir sfc albedo for diffuse radiation (rsirdf) fraction m*ndim
955 !
956 !----- Output parameters
957 !
958 ! all-sky flux (downward minus upward) (flx) fraction m*ndim*(np+1)
959 ! clear-sky flux (downward minus upward) (flc) fraction m*ndim*(np+1)
960 ! all-sky direct downward uv (0.175-0.4 micron)
961 ! flux at the surface (fdiruv) fraction m*ndim
962 ! all-sky diffuse downward uv flux at
963 ! the surface (fdifuv) fraction m*ndim
964 ! all-sky direct downward par (0.4-0.7 micron)
965 ! flux at the surface (fdirpar) fraction m*ndim
966 ! all-sky diffuse downward par flux at
967 ! the surface (fdifpar) fraction m*ndim
968 ! all-sky direct downward ir (0.7-10 micron)
969 ! flux at the surface (fdirir) fraction m*ndim
970 ! all-sky diffuse downward ir flux at
971 ! the surface (fdifir) fraction m*ndim
972 !
973 !----- Notes:
974 !
975 ! (1) The unit of "flux" is fraction of the incoming solar radiation
976 ! at the top of the atmosphere. Therefore, fluxes should
977 ! be equal to "flux" multiplied by the extra-terrestrial solar
978 ! flux and the cosine of solar zenith angle.
979 ! (2) pl(i,j,1) is the pressure at the top of the model, and
980 ! pl(i,j,np+1) is the surface pressure.
981 ! (3) the pressure levels ict and icb correspond approximately
982 ! to 400 and 700 mb.
983 ! (4) if overcast='true', the clear-sky flux, flc, is not computed.
984 !
985 !**************************************************************************
986 implicit none
987 !**************************************************************************
988
989 !-----input parameters
990
991 integer m,n,ndim,np
992 integer ict(m,ndim),icb(m,ndim)
993 real pl(m,ndim,np+1),ta(m,ndim,np),wa(m,ndim,np),oa(m,ndim,np)
994 real cwc(m,ndim,np,2),taucld(m,ndim,np,2),reff(m,ndim,np,2), &
995 fcld(m,ndim,np)
996 real taual(m,ndim,np,11),ssaal(m,ndim,np,11),asyal(m,ndim,np,11)
997 real cosz(m,ndim),rsuvbm(m,ndim),rsuvdf(m,ndim), &
998 rsirbm(m,ndim),rsirdf(m,ndim)
999 logical overcast,cldwater
1000
1001 !-----output parameters
1002
1003 real flx(m,ndim,np+1),flc(m,ndim,np+1)
1004 real flxu(m,ndim,np+1),flxd(m,ndim,np+1)
1005 real fdiruv (m,ndim),fdifuv (m,ndim)
1006 real fdirpar(m,ndim),fdifpar(m,ndim)
1007 real fdirir (m,ndim),fdifir (m,ndim)
1008
1009 !-----temporary array
1010
1011 integer i,j,k
1012 real cwp(m,n,np,2)
1013 real dp(m,n,np),wh(m,n,np),oh(m,n,np),scal(m,n,np)
1014 real swh(m,n,np+1),so2(m,n,np+1),df(m,n,np+1)
1015 real sdf(m,n),sclr(m,n),csm(m,n),x
1016
1017 do j= 1, n
1018 do i= 1, m
1019 if (pl(i,j,1) .eq. 0.0) then
1020 pl(i,j,1)=1.0e-4
1021 endif
1022 enddo
1023 enddo
1024
1025 do j= 1, n
1026 do i= 1, m
1027
1028 swh(i,j,1)=0.
1029 so2(i,j,1)=0.
1030
1031 !-----csm is the effective secant of the solar zenith angle
1032 ! see equation (12) of Lacis and Hansen (1974, JAS)
1033
1034 csm(i,j)=35./sqrt(1224.*cosz(i,j)*cosz(i,j)+1.)
1035
1036 enddo
1037 enddo
1038
1039 do k= 1, np
1040 do j= 1, n
1041 do i= 1, m
1042
1043 !-----compute layer thickness and pressure-scaling function.
1044 ! indices for the surface level and surface layer
1045 ! are np+1 and np, respectively.
1046
1047 dp(i,j,k)=pl(i,j,k+1)-pl(i,j,k)
1048 scal(i,j,k)=dp(i,j,k)*(.5*(pl(i,j,k)+pl(i,j,k+1))/300.)**.8
1049
1050 !-----compute scaled water vapor amount, unit is g/cm**2
1051 ! note: the sign prior to the constant 0.00135 was incorrectly
1052 ! set to negative in the previous version
1053
1054 wh(i,j,k)=1.02*wa(i,j,k)*scal(i,j,k)* &
1055 (1.+0.00135*(ta(i,j,k)-240.)) +1.e-11
1056 swh(i,j,k+1)=swh(i,j,k)+wh(i,j,k)
1057
1058 !-----compute ozone amount, unit is (cm-atm)stp
1059 ! the number 466.7 is a conversion factor from g/cm**2 to (cm-atm)stp
1060
1061 oh(i,j,k)=1.02*oa(i,j,k)*dp(i,j,k)*466.7 +1.e-11
1062
1063 !-----compute layer cloud water amount (gm/m**2)
1064 ! the index is 1 for ice crystals and 2 for liquid drops
1065
1066 cwp(i,j,k,1)=1.02*10000.*cwc(i,j,k,1)*dp(i,j,k)
1067 cwp(i,j,k,2)=1.02*10000.*cwc(i,j,k,2)*dp(i,j,k)
1068
1069 enddo
1070 enddo
1071 enddo
1072
1073 !-----initialize fluxes for all-sky (flx), clear-sky (flc), and
1074 ! flux reduction (df)
1075
1076 do k=1, np+1
1077 do j=1, n
1078 do i=1, m
1079 flx(i,j,k)=0.
1080 flc(i,j,k)=0.
1081 flxu(i,j,k)=0.
1082 flxd(i,j,k)=0.
1083 df(i,j,k)=0.
1084 enddo
1085 enddo
1086 enddo
1087
1088 !-----compute solar uv and par fluxes
1089
1090 call soluv (m,n,ndim,np,oh,dp,overcast,cldwater, &
1091 cwp,taucld,reff,ict,icb,fcld,cosz, &
1092 taual,ssaal,asyal,csm,rsuvbm,rsuvdf, &
1093 flx,flc,flxu,flxd,fdiruv,fdifuv,fdirpar,fdifpar)
1094
1095 !-----compute and update solar ir fluxes
1096
1097 call solir (m,n,ndim,np,wh,overcast,cldwater, &
1098 cwp,taucld,reff,ict,icb,fcld,cosz, &
1099 taual,ssaal,asyal,csm,rsirbm,rsirdf, &
1100 flx,flc,flxu,flxd,fdirir,fdifir)
1101
1102 !-----compute scaled o2 amount, unit is (cm-atm)stp.
1103
1104 do k= 1, np
1105 do j= 1, n
1106 do i= 1, m
1107 so2(i,j,k+1)=so2(i,j,k)+165.22*scal(i,j,k)
1108 enddo
1109 enddo
1110 enddo
1111
1112 !-----compute flux reduction due to oxygen following
1113 ! chou (J. climate, 1990). The fraction 0.0287 is the
1114 ! extraterrestrial solar flux in the o2 bands.
1115
1116 do k= 2, np+1
1117 do j= 1, n
1118 do i= 1, m
1119 x=so2(i,j,k)*csm(i,j)
1120 df(i,j,k)=df(i,j,k)+0.0287*(1.-exp(-0.00027*sqrt(x)))
1121 enddo
1122 enddo
1123 enddo
1124
1125 !-----compute scaled co2 amounts. unit is (cm-atm)stp.
1126
1127 do k= 1, np
1128 do j= 1, n
1129 do i= 1, m
1130 so2(i,j,k+1)=so2(i,j,k)+co2*789.*scal(i,j,k)+1.e-11
1131 enddo
1132 enddo
1133 enddo
1134
1135 !-----compute and update flux reduction due to co2 following
1136 ! chou (J. Climate, 1990)
1137
1138 call flxco2(m,n,np,so2,swh,csm,df)
1139
1140 !-----adjust for the effect of o2 cnd co2 on clear-sky fluxes.
1141
1142 do k= 2, np+1
1143 do j= 1, n
1144 do i= 1, m
1145 flc(i,j,k)=flc(i,j,k)-df(i,j,k)
1146 enddo
1147 enddo
1148 enddo
1149
1150 !-----adjust for the all-sky fluxes due to o2 and co2. It is
1151 ! assumed that o2 and co2 have no effects on solar radiation
1152 ! below clouds.
1153
1154 do j=1,n
1155 do i=1,m
1156 sdf(i,j)=0.0
1157 sclr(i,j)=1.0
1158 enddo
1159 enddo
1160
1161 do k=1,np
1162 do j=1,n
1163 do i=1,m
1164
1165 !-----sclr is the fraction of clear sky.
1166 ! sdf is the flux reduction below clouds.
1167
1168 if(fcld(i,j,k).gt.0.01) then
1169 sdf(i,j)=sdf(i,j)+df(i,j,k)*sclr(i,j)*fcld(i,j,k)
1170 sclr(i,j)=sclr(i,j)*(1.-fcld(i,j,k))
1171 endif
1172 flx(i,j,k+1)=flx(i,j,k+1)-sdf(i,j)-df(i,j,k+1)*sclr(i,j)
1173 flxu(i,j,k+1)=flxu(i,j,k+1)-sdf(i,j)-df(i,j,k+1)*sclr(i,j)
1174 flxd(i,j,k+1)=flxd(i,j,k+1)-sdf(i,j)-df(i,j,k+1)*sclr(i,j) ! SG: same as flux????
1175
1176 enddo
1177 enddo
1178 enddo
1179
1180 !-----adjustment for the direct downward ir flux.
1181
1182 do j= 1, n
1183 do i= 1, m
1184 flc(i,j,np+1)=flc(i,j,np+1)+df(i,j,np+1)*rsirbm(i,j)
1185 flx(i,j,np+1)=flx(i,j,np+1)+(sdf(i,j)+ &
1186 df(i,j,np+1)*sclr(i,j))*rsirbm(i,j)
1187 flxu(i,j,np+1)=flxu(i,j,np+1)+(sdf(i,j)+ &
1188 df(i,j,np+1)*sclr(i,j))*rsirbm(i,j)
1189 flxd(i,j,np+1)=flxd(i,j,np+1)+(sdf(i,j)+ &
1190 df(i,j,np+1)*sclr(i,j))*rsirbm(i,j)
1191 fdirir(i,j)=fdirir(i,j)-(sdf(i,j)+df(i,j,np+1)*sclr(i,j))
1192 enddo
1193 enddo
1194
1195 end subroutine sorad
1196
1197 !************************************************************************
1198
1199 subroutine soluv (m,n,ndim,np,oh,dp,overcast,cldwater, &
1200 cwp,taucld,reff,ict,icb,fcld,cosz, &
1201 taual,ssaal,asyal,csm,rsuvbm,rsuvdf, &
1202 flx,flc,flxu,flxd,fdiruv,fdifuv,fdirpar,fdifpar)
1203
1204 !************************************************************************
1205 ! compute solar fluxes in the uv+par region. the spectrum is
1206 ! grouped into 8 bands:
1207 !
1208 ! Band Micrometer
1209 !
1210 ! UV-C 1. .175 - .225
1211 ! 2. .225 - .245
1212 ! .260 - .280
1213 ! 3. .245 - .260
1214 !
1215 ! UV-B 4. .280 - .295
1216 ! 5. .295 - .310
1217 ! 6. .310 - .320
1218 !
1219 ! UV-A 7. .320 - .400
1220 !
1221 ! PAR 8. .400 - .700
1222 !
1223 !----- Input parameters: units size
1224 !
1225 ! number of soundings in zonal direction (m) n/d 1
1226 ! number of soundings in meridional direction (n) n/d 1
1227 ! maximum number of soundings in n/d 1
1228 ! meridional direction (ndim)
1229 ! number of atmospheric layers (np) n/d 1
1230 ! layer ozone content (oh) (cm-atm)stp m*n*np
1231 ! layer pressure thickness (dp) mb m*n*np
1232 ! option for scaling cloud optical thickness n/d 1
1233 ! overcast="true" if scaling is NOT required
1234 ! overcast="fasle" if scaling is required
1235 ! input option for cloud optical thickness n/d 1
1236 ! cldwater="true" if taucld is provided
1237 ! cldwater="false" if cwp is provided
1238 ! cloud water amount (cwp) gm/m**2 m*n*np*2
1239 ! index 1 for ice particles
1240 ! index 2 for liquid drops
1241 ! cloud optical thickness (taucld) n/d m*ndim*np*2
1242 ! index 1 for ice paticles
1243 ! index 2 for liquid particles
1244 ! effective cloud-particle size (reff) micrometer m*ndim*np*2
1245 ! index 1 for ice paticles
1246 ! index 2 for liquid particles
1247 ! level indiex separating high and n/d m*n
1248 ! middle clouds (ict)
1249 ! level indiex separating middle and n/d m*n
1250 ! low clouds (icb)
1251 ! cloud amount (fcld) fraction m*ndim*np
1252 ! cosine of solar zenith angle (cosz) n/d m*ndim
1253 ! aerosol optical thickness (taual) n/d m*ndim*np*11
1254 ! aerosol single-scattering albedo (ssaal) n/d m*ndim*np*11
1255 ! aerosol asymmetry factor (asyal) n/d m*ndim*np*11
1256 ! cosecant of the solar zenith angle (csm) n/d m*n
1257 ! uv+par surface albedo for beam fraction m*ndim
1258 ! radiation (rsuvbm)
1259 ! uv+par surface albedo for diffuse fraction m*ndim
1260 ! radiation (rsuvdf)
1261 !
1262 !---- temporary array
1263 !
1264 ! scaled cloud optical thickness n/d m*n*np
1265 ! for beam radiation (tauclb)
1266 ! scaled cloud optical thickness n/d m*n*np
1267 ! for diffuse radiation (tauclf)
1268 !
1269 !----- output (updated) parameters:
1270 !
1271 ! all-sky net downward flux (flx) fraction m*ndim*(np+1)
1272 ! clear-sky net downward flux (flc) fraction m*ndim*(np+1)
1273 ! all-sky direct downward uv flux at
1274 ! the surface (fdiruv) fraction m*ndim
1275 ! all-sky diffuse downward uv flux at
1276 ! the surface (fdifuv) fraction m*ndim
1277 ! all-sky direct downward par flux at
1278 ! the surface (fdirpar) fraction m*ndim
1279 ! all-sky diffuse downward par flux at
1280 ! the surface (fdifpar) fraction m*ndim
1281 !
1282 !***********************************************************************
1283 implicit none
1284 !***********************************************************************
1285
1286 !-----input parameters
1287
1288 integer m,n,ndim,np
1289 integer ict(m,ndim),icb(m,ndim)
1290 real taucld(m,ndim,np,2),reff(m,ndim,np,2),fcld(m,ndim,np)
1291 real cc(m,n,3),cosz(m,ndim)
1292 real cwp(m,n,np,2),oh(m,n,np),dp(m,n,np)
1293 real taual(m,ndim,np,11),ssaal(m,ndim,np,11),asyal(m,ndim,np,11)
1294 real rsuvbm(m,ndim),rsuvdf(m,ndim),csm(m,n)
1295 logical overcast,cldwater
1296
1297 !-----output (updated) parameter
1298
1299 real flx(m,ndim,np+1),flc(m,ndim,np+1)
1300 real flxu(m,ndim,np+1),flxd(m,ndim,np+1)
1301 real fdiruv (m,ndim),fdifuv (m,ndim)
1302 real fdirpar(m,ndim),fdifpar(m,ndim)
1303
1304 !-----static parameters
1305
1306 integer nband
1307 parameter (nband=8)
1308 real hk(nband),xk(nband),ry(nband)
1309 real aig(3),awg(3)
1310
1311 !-----temporary array
1312
1313 integer i,j,k,ib
1314 real tauclb(m,n,np),tauclf(m,n,np),asycl(m,n,np)
1315 real taurs,tauoz,tausto,ssatau,asysto,tauto,ssato,asyto
1316 real taux,reff1,reff2,g1,g2
1317 real td(m,n,np+1,2),rr(m,n,np+1,2),tt(m,n,np+1,2), &
1318 rs(m,n,np+1,2),ts(m,n,np+1,2)
1319 real fall(m,n,np+1),fclr(m,n,np+1),fsdir(m,n),fsdif(m,n)
1320 real fallu(m,n,np+1),falld(m,n,np+1)
1321 real asyclt(m,n)
1322 real rr1t(m,n),tt1t(m,n),td1t(m,n),rs1t(m,n),ts1t(m,n)
1323 real rr2t(m,n),tt2t(m,n),td2t(m,n),rs2t(m,n),ts2t(m,n)
1324
1325 !-----hk is the fractional extra-terrestrial solar flux in each
1326 ! of the 8 bands. the sum of hk is 0.47074.
1327
1328 data hk/.00057, .00367, .00083, .00417, &
1329 .00600, .00556, .05913, .39081/
1330
1331 !-----xk is the ozone absorption coefficient. unit: /(cm-atm)stp
1332
1333 data xk /30.47, 187.2, 301.9, 42.83, &
1334 7.09, 1.25, 0.0345, 0.0539/
1335
1336 !-----ry is the extinction coefficient for Rayleigh scattering.
1337 ! unit: /mb.
1338
1339 data ry /.00604, .00170, .00222, .00132, &
1340 .00107, .00091, .00055, .00012/
1341
1342 !-----coefficients for computing the asymmetry factor of ice clouds
1343 ! from asycl=aig(*,1)+aig(*,2)*reff+aig(*,3)*reff**2, independent
1344 ! of spectral band.
1345
1346 data aig/.74625000,.00105410,-.00000264/
1347
1348 !-----coefficients for computing the asymmetry factor of liquid
1349 ! clouds from asycl=awg(*,1)+awg(*,2)*reff+awg(*,3)*reff**2,
1350 ! independent of spectral band.
1351
1352 data awg/.82562000,.00529000,-.00014866/
1353
1354 !-----initialize fdiruv, fdifuv, surface reflectances and transmittances.
1355 ! cc is the maximum cloud cover in each of the three cloud groups.
1356
1357 do j= 1, n
1358 do i= 1, m
1359 fdiruv(i,j)=0.0
1360 fdifuv(i,j)=0.0
1361 rr(i,j,np+1,1)=rsuvbm(i,j)
1362 rr(i,j,np+1,2)=rsuvbm(i,j)
1363 rs(i,j,np+1,1)=rsuvdf(i,j)
1364 rs(i,j,np+1,2)=rsuvdf(i,j)
1365 td(i,j,np+1,1)=0.0
1366 td(i,j,np+1,2)=0.0
1367 tt(i,j,np+1,1)=0.0
1368 tt(i,j,np+1,2)=0.0
1369 ts(i,j,np+1,1)=0.0
1370 ts(i,j,np+1,2)=0.0
1371 cc(i,j,1)=0.0
1372 cc(i,j,2)=0.0
1373 cc(i,j,3)=0.0
1374 enddo
1375 enddo
1376
1377
1378 !-----compute cloud optical thickness
1379
1380 if (cldwater) then
1381
1382 do k= 1, np
1383 do j= 1, n
1384 do i= 1, m
1385 taucld(i,j,k,1)=cwp(i,j,k,1)*( 3.33e-4+2.52/reff(i,j,k,1))
1386 taucld(i,j,k,2)=cwp(i,j,k,2)*(-6.59e-3+1.65/reff(i,j,k,2))
1387 enddo
1388 enddo
1389 enddo
1390
1391 endif
1392
1393 !-----options for scaling cloud optical thickness
1394
1395 if (overcast) then
1396
1397 do k= 1, np
1398 do j= 1, n
1399 do i= 1, m
1400 tauclb(i,j,k)=taucld(i,j,k,1)+taucld(i,j,k,2)
1401 tauclf(i,j,k)=tauclb(i,j,k)
1402 enddo
1403 enddo
1404 enddo
1405
1406 do k= 1, 3
1407 do j= 1, n
1408 do i= 1, m
1409 cc(i,j,k)=1.0
1410 enddo
1411 enddo
1412 enddo
1413
1414 else
1415
1416 !-----scale cloud optical thickness in each layer from taucld (with
1417 ! cloud amount fcld) to tauclb and tauclf (with cloud amount cc).
1418 ! tauclb is the scaled optical thickness for beam radiation and
1419 ! tauclf is for diffuse radiation.
1420
1421 call cldscale(m,n,ndim,np,cosz,fcld,taucld,ict,icb, &
1422 cc,tauclb,tauclf)
1423
1424 endif
1425
1426 !-----compute cloud asymmetry factor for a mixture of
1427 ! liquid and ice particles. unit of reff is micrometers.
1428
1429 do k= 1, np
1430
1431 do j= 1, n
1432 do i= 1, m
1433
1434 asyclt(i,j)=1.0
1435
1436 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1437 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1438
1439 reff1=min(reff(i,j,k,1),130.)
1440 reff2=min(reff(i,j,k,2),20.0)
1441
1442 g1=(aig(1)+(aig(2)+aig(3)*reff1)*reff1)*taucld(i,j,k,1)
1443 g2=(awg(1)+(awg(2)+awg(3)*reff2)*reff2)*taucld(i,j,k,2)
1444 asyclt(i,j)=(g1+g2)/taux
1445
1446 endif
1447
1448 enddo
1449 enddo
1450
1451 do j=1,n
1452 do i=1,m
1453 asycl(i,j,k)=asyclt(i,j)
1454 enddo
1455 enddo
1456
1457 enddo
1458
1459 !-----integration over spectral bands
1460
1461 do 100 ib=1,nband
1462
1463 do 300 k= 1, np
1464
1465 do j= 1, n
1466 do i= 1, m
1467
1468 !-----compute ozone and rayleigh optical thicknesses
1469
1470 taurs=ry(ib)*dp(i,j,k)
1471 tauoz=xk(ib)*oh(i,j,k)
1472
1473 !-----compute clear-sky optical thickness, single scattering albedo,
1474 ! and asymmetry factor
1475
1476 tausto=taurs+tauoz+taual(i,j,k,ib)+1.0e-8
1477 ssatau=ssaal(i,j,k,ib)*taual(i,j,k,ib)+taurs
1478 asysto=asyal(i,j,k,ib)*ssaal(i,j,k,ib)*taual(i,j,k,ib)
1479
1480 tauto=tausto
1481 ssato=ssatau/tauto+1.0e-8
1482 ssato=min(ssato,0.999999)
1483 asyto=asysto/(ssato*tauto)
1484
1485 !-----compute reflectance and transmittance for cloudless layers
1486
1487 !- for direct incident radiation
1488
1489 call deledd (tauto,ssato,asyto,csm(i,j), &
1490 rr1t(i,j),tt1t(i,j),td1t(i,j))
1491
1492 !- for diffuse incident radiation
1493
1494 call sagpol (tauto,ssato,asyto,rs1t(i,j),ts1t(i,j))
1495
1496 !-----compute reflectance and transmittance for cloud layers
1497
1498 if (tauclb(i,j,k).lt.0.01 .or. fcld(i,j,k).lt.0.01) then
1499
1500 rr2t(i,j)=rr1t(i,j)
1501 tt2t(i,j)=tt1t(i,j)
1502 td2t(i,j)=td1t(i,j)
1503 rs2t(i,j)=rs1t(i,j)
1504 ts2t(i,j)=ts1t(i,j)
1505
1506 else
1507
1508 !-- for direct incident radiation
1509
1510 tauto=tausto+tauclb(i,j,k)
1511 ssato=(ssatau+tauclb(i,j,k))/tauto+1.0e-8
1512 ssato=min(ssato,0.999999)
1513 asyto=(asysto+asycl(i,j,k)*tauclb(i,j,k))/(ssato*tauto)
1514
1515 call deledd (tauto,ssato,asyto,csm(i,j), &
1516 rr2t(i,j),tt2t(i,j),td2t(i,j))
1517
1518 !-- for diffuse incident radiation
1519
1520 tauto=tausto+tauclf(i,j,k)
1521 ssato=(ssatau+tauclf(i,j,k))/tauto+1.0e-8
1522 ssato=min(ssato,0.999999)
1523 asyto=(asysto+asycl(i,j,k)*tauclf(i,j,k))/(ssato*tauto)
1524
1525 call sagpol (tauto,ssato,asyto,rs2t(i,j),ts2t(i,j))
1526
1527 endif
1528
1529 enddo
1530 enddo
1531
1532 do j=1,n
1533 do i=1,m
1534 rr(i,j,k,1)=rr1t(i,j)
1535 enddo
1536 enddo
1537 do j=1,n
1538 do i=1,m
1539 tt(i,j,k,1)=tt1t(i,j)
1540 enddo
1541 enddo
1542 do j=1,n
1543 do i=1,m
1544 td(i,j,k,1)=td1t(i,j)
1545 enddo
1546 enddo
1547 do j=1,n
1548 do i=1,m
1549 rs(i,j,k,1)=rs1t(i,j)
1550 enddo
1551 enddo
1552 do j=1,n
1553 do i=1,m
1554 ts(i,j,k,1)=ts1t(i,j)
1555 enddo
1556 enddo
1557
1558 do j=1,n
1559 do i=1,m
1560 rr(i,j,k,2)=rr2t(i,j)
1561 enddo
1562 enddo
1563 do j=1,n
1564 do i=1,m
1565 tt(i,j,k,2)=tt2t(i,j)
1566 enddo
1567 enddo
1568 do j=1,n
1569 do i=1,m
1570 td(i,j,k,2)=td2t(i,j)
1571 enddo
1572 enddo
1573 do j=1,n
1574 do i=1,m
1575 rs(i,j,k,2)=rs2t(i,j)
1576 enddo
1577 enddo
1578 do j=1,n
1579 do i=1,m
1580 ts(i,j,k,2)=ts2t(i,j)
1581 enddo
1582 enddo
1583
1584 300 continue
1585
1586 !-----flux calculations
1587
1588 call cldflx (m,n,np,ict,icb,overcast,cc,rr,tt,td,rs,ts, &
1589 fclr,fall,fallu,falld,fsdir,fsdif)
1590
1591 do k= 1, np+1
1592 do j= 1, n
1593 do i= 1, m
1594 flx(i,j,k)=flx(i,j,k)+fall(i,j,k)*hk(ib)
1595 flxu(i,j,k)=flxu(i,j,k)+fallu(i,j,k)*hk(ib)
1596 flxd(i,j,k)=flxd(i,j,k)+falld(i,j,k)*hk(ib)
1597 enddo
1598 enddo
1599 do j= 1, n
1600 do i= 1, m
1601 flc(i,j,k)=flc(i,j,k)+fclr(i,j,k)*hk(ib)
1602 enddo
1603 enddo
1604 enddo
1605
1606 !-----compute downward surface fluxes in the UV and par regions
1607
1608 if(ib.lt.8) then
1609 do j=1,n
1610 do i=1,m
1611 fdiruv(i,j)=fdiruv(i,j)+fsdir(i,j)*hk(ib)
1612 fdifuv(i,j)=fdifuv(i,j)+fsdif(i,j)*hk(ib)
1613 enddo
1614 enddo
1615 else
1616 do j=1,n
1617 do i=1,m
1618 fdirpar(i,j)=fsdir(i,j)*hk(ib)
1619 fdifpar(i,j)=fsdif(i,j)*hk(ib)
1620 enddo
1621 enddo
1622 endif
1623
1624 100 continue
1625
1626 end subroutine soluv
1627
1628 !************************************************************************
1629
1630 subroutine solir (m,n,ndim,np,wh,overcast,cldwater, &
1631 cwp,taucld,reff,ict,icb,fcld,cosz, &
1632 taual,ssaal,asyal,csm,rsirbm,rsirdf, &
1633 flx,flc,flxu,flxd,fdirir,fdifir)
1634
1635 !************************************************************************
1636 ! compute solar flux in the infrared region. The spectrum is divided
1637 ! into three bands:
1638 !
1639 ! band wavenumber(/cm) wavelength (micron)
1640 ! 1( 9) 14300-8200 0.70-1.22
1641 ! 2(10) 8200-4400 1.22-2.27
1642 ! 3(11) 4400-1000 2.27-10.0
1643 !
1644 !----- Input parameters: units size
1645 !
1646 ! number of soundings in zonal direction (m) n/d 1
1647 ! number of soundings in meridional direction (n) n/d 1
1648 ! maximum number of soundings in n/d 1
1649 ! meridional direction (ndim)
1650 ! number of atmospheric layers (np) n/d 1
1651 ! layer scaled-water vapor content (wh) gm/cm^2 m*n*np
1652 ! option for scaling cloud optical thickness n/d 1
1653 ! overcast="true" if scaling is NOT required
1654 ! overcast="fasle" if scaling is required
1655 ! input option for cloud optical thickness n/d 1
1656 ! cldwater="true" if taucld is provided
1657 ! cldwater="false" if cwp is provided
1658 ! cloud water concentration (cwp) gm/m**2 m*n*np*2
1659 ! index 1 for ice particles
1660 ! index 2 for liquid drops
1661 ! cloud optical thickness (taucld) n/d m*ndim*np*2
1662 ! index 1 for ice paticles
1663 ! effective cloud-particle size (reff) micrometer m*ndim*np*2
1664 ! index 1 for ice paticles
1665 ! index 2 for liquid particles
1666 ! level index separating high and n/d m*n
1667 ! middle clouds (ict)
1668 ! level index separating middle and n/d m*n
1669 ! low clouds (icb)
1670 ! cloud amount (fcld) fraction m*ndim*np
1671 ! aerosol optical thickness (taual) n/d m*ndim*np*11
1672 ! aerosol single-scattering albedo (ssaal) n/d m*ndim*np*11
1673 ! aerosol asymmetry factor (asyal) n/d m*ndim*np*11
1674 ! cosecant of the solar zenith angle (csm) n/d m*n
1675 ! near ir surface albedo for beam fraction m*ndim
1676 ! radiation (rsirbm)
1677 ! near ir surface albedo for diffuse fraction m*ndim
1678 ! radiation (rsirdf)
1679 !
1680 !---- temporary array
1681 !
1682 ! scaled cloud optical thickness n/d m*n*np
1683 ! for beam radiation (tauclb)
1684 ! scaled cloud optical thickness n/d m*n*np
1685 ! for diffuse radiation (tauclf)
1686 !
1687 !----- output (updated) parameters:
1688 !
1689 ! all-sky flux (downward-upward) (flx) fraction m*ndim*(np+1)
1690 ! clear-sky flux (downward-upward) (flc) fraction m*ndim*(np+1)
1691 ! all-sky direct downward ir flux at
1692 ! the surface (fdirir) fraction m*ndim
1693 ! all-sky diffuse downward ir flux at
1694 ! the surface (fdifir) fraction m*ndim
1695 !
1696 !**********************************************************************
1697 implicit none
1698 !**********************************************************************
1699
1700 !-----input parameters
1701
1702 integer m,n,ndim,np
1703 integer ict(m,ndim),icb(m,ndim)
1704 real cwp(m,n,np,2),taucld(m,ndim,np,2),reff(m,ndim,np,2)
1705 real fcld(m,ndim,np),cc(m,n,3),cosz(m,ndim)
1706 real rsirbm(m,ndim),rsirdf(m,ndim)
1707 real taual(m,ndim,np,11),ssaal(m,ndim,np,11),asyal(m,ndim,np,11)
1708 real wh(m,n,np),csm(m,n)
1709 logical overcast,cldwater
1710
1711 !-----output (updated) parameters
1712
1713 real flx(m,ndim,np+1),flc(m,ndim,np+1)
1714 real flxu(m,ndim,np+1),flxd(m,ndim,np+1)
1715 real fdirir(m,ndim),fdifir(m,ndim)
1716
1717 !-----static parameters
1718
1719 integer nk,nband
1720 parameter (nk=10,nband=3)
1721 real xk(nk),hk(nband,nk),aib(nband,2),awb(nband,2)
1722 real aia(nband,3),awa(nband,3),aig(nband,3),awg(nband,3)
1723
1724 !-----temporary array
1725
1726 integer ib,iv,ik,i,j,k
1727 real tauclb(m,n,np),tauclf(m,n,np)
1728 real ssacl(m,n,np),asycl(m,n,np)
1729 real rr(m,n,np+1,2),tt(m,n,np+1,2),td(m,n,np+1,2), &
1730 rs(m,n,np+1,2),ts(m,n,np+1,2)
1731 real fall(m,n,np+1),fclr(m,n,np+1)
1732 real fallu(m,n,np+1),falld(m,n,np+1)
1733 real fsdir(m,n),fsdif(m,n)
1734
1735 real tauwv,tausto,ssatau,asysto,tauto,ssato,asyto
1736 real taux,reff1,reff2,w1,w2,g1,g2
1737 real ssaclt(m,n),asyclt(m,n)
1738 real rr1t(m,n),tt1t(m,n),td1t(m,n),rs1t(m,n),ts1t(m,n)
1739 real rr2t(m,n),tt2t(m,n),td2t(m,n),rs2t(m,n),ts2t(m,n)
1740
1741 !-----water vapor absorption coefficient for 10 k-intervals.
1742 ! unit: cm^2/gm
1743
1744 data xk/ &
1745 0.0010, 0.0133, 0.0422, 0.1334, 0.4217, &
1746 1.334, 5.623, 31.62, 177.8, 1000.0/
1747
1748 !-----water vapor k-distribution function,
1749 ! the sum of hk is 0.52926. unit: fraction
1750
1751 data hk/ &
1752 .20673,.08236,.01074, .03497,.01157,.00360, &
1753 .03011,.01133,.00411, .02260,.01143,.00421, &
1754 .01336,.01240,.00389, .00696,.01258,.00326, &
1755 .00441,.01381,.00499, .00115,.00650,.00465, &
1756 .00026,.00244,.00245, .00000,.00094,.00145/
1757
1758 !-----coefficients for computing the extinction coefficient of
1759 ! ice clouds from b=aib(*,1)+aib(*,2)/reff
1760
1761 data aib/ &
1762 .000333, .000333, .000333, &
1763 2.52, 2.52, 2.52/
1764
1765 !-----coefficients for computing the extinction coefficient of
1766 ! water clouds from b=awb(*,1)+awb(*,2)/reff
1767
1768 data awb/ &
1769 -0.0101, -0.0166, -0.0339, &
1770 1.72, 1.85, 2.16/
1771
1772
1773 !-----coefficients for computing the single scattering albedo of
1774 ! ice clouds from ssa=1-(aia(*,1)+aia(*,2)*reff+aia(*,3)*reff**2)
1775
1776 data aia/ &
1777 -.00000260, .00215346, .08938331, &
1778 .00000746, .00073709, .00299387, &
1779 .00000000,-.00000134,-.00001038/
1780
1781 !-----coefficients for computing the single scattering albedo of
1782 ! liquid clouds from ssa=1-(awa(*,1)+awa(*,2)*reff+awa(*,3)*reff**2)
1783
1784 data awa/ &
1785 .00000007,-.00019934, .01209318, &
1786 .00000845, .00088757, .01784739, &
1787 -.00000004,-.00000650,-.00036910/
1788
1789 !-----coefficients for computing the asymmetry factor of ice clouds
1790 ! from asycl=aig(*,1)+aig(*,2)*reff+aig(*,3)*reff**2
1791
1792 data aig/ &
1793 .74935228, .76098937, .84090400, &
1794 .00119715, .00141864, .00126222, &
1795 -.00000367,-.00000396,-.00000385/
1796
1797 !-----coefficients for computing the asymmetry factor of liquid clouds
1798 ! from asycl=awg(*,1)+awg(*,2)*reff+awg(*,3)*reff**2
1799
1800 data awg/ &
1801 .79375035, .74513197, .83530748, &
1802 .00832441, .01370071, .00257181, &
1803 -.00023263,-.00038203, .00005519/
1804
1805 !-----initialize surface fluxes, reflectances, and transmittances.
1806 ! cc is the maximum cloud cover in each of the three cloud groups.
1807
1808 do j= 1, n
1809 do i= 1, m
1810 fdirir(i,j)=0.0
1811 fdifir(i,j)=0.0
1812 rr(i,j,np+1,1)=rsirbm(i,j)
1813 rr(i,j,np+1,2)=rsirbm(i,j)
1814 rs(i,j,np+1,1)=rsirdf(i,j)
1815 rs(i,j,np+1,2)=rsirdf(i,j)
1816 td(i,j,np+1,1)=0.0
1817 td(i,j,np+1,2)=0.0
1818 tt(i,j,np+1,1)=0.0
1819 tt(i,j,np+1,2)=0.0
1820 ts(i,j,np+1,1)=0.0
1821 ts(i,j,np+1,2)=0.0
1822 cc(i,j,1)=0.0
1823 cc(i,j,2)=0.0
1824 cc(i,j,3)=0.0
1825 enddo
1826 enddo
1827
1828 !-----integration over spectral bands
1829
1830 do 100 ib=1,nband
1831
1832 iv=ib+8
1833
1834 !-----compute cloud optical thickness
1835
1836 if (cldwater) then
1837
1838 do k= 1, np
1839 do j= 1, n
1840 do i= 1, m
1841 taucld(i,j,k,1)=cwp(i,j,k,1)*(aib(ib,1) &
1842 +aib(ib,2)/reff(i,j,k,1))
1843 taucld(i,j,k,2)=cwp(i,j,k,2)*(awb(ib,1) &
1844 +awb(ib,2)/reff(i,j,k,2))
1845 enddo
1846 enddo
1847 enddo
1848
1849 endif
1850
1851 !-----options for scaling cloud optical thickness
1852
1853 if (overcast) then
1854
1855 do k= 1, np
1856 do j= 1, n
1857 do i= 1, m
1858 tauclb(i,j,k)=taucld(i,j,k,1)+taucld(i,j,k,2)
1859 tauclf(i,j,k)=tauclb(i,j,k)
1860 enddo
1861 enddo
1862 enddo
1863
1864 do k= 1, 3
1865 do j= 1, n
1866 do i= 1, m
1867 cc(i,j,k)=1.0
1868 enddo
1869 enddo
1870 enddo
1871
1872 else
1873
1874 !-----scale cloud optical thickness in each layer from taucld (with
1875 ! cloud amount fcld) to tauclb and tauclf (with cloud amount cc).
1876 ! tauclb is the scaled optical thickness for beam radiation and
1877 ! tauclf is for diffuse radiation.
1878
1879 call cldscale(m,n,ndim,np,cosz,fcld,taucld,ict,icb, &
1880 cc,tauclb,tauclf)
1881
1882 endif
1883
1884 !-----compute cloud single scattering albedo and asymmetry factor
1885 ! for a mixture of ice and liquid particles.
1886
1887 do k= 1, np
1888
1889 do j= 1, n
1890 do i= 1, m
1891
1892 ssaclt(i,j)=1.0
1893 asyclt(i,j)=1.0
1894
1895 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1896 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1897
1898 reff1=min(reff(i,j,k,1),130.)
1899 reff2=min(reff(i,j,k,2),20.0)
1900
1901 w1=(1.-(aia(ib,1)+(aia(ib,2)+ &
1902 aia(ib,3)*reff1)*reff1))*taucld(i,j,k,1)
1903 w2=(1.-(awa(ib,1)+(awa(ib,2)+ &
1904 awa(ib,3)*reff2)*reff2))*taucld(i,j,k,2)
1905 ssaclt(i,j)=(w1+w2)/taux
1906
1907 g1=(aig(ib,1)+(aig(ib,2)+aig(ib,3)*reff1)*reff1)*w1
1908 g2=(awg(ib,1)+(awg(ib,2)+awg(ib,3)*reff2)*reff2)*w2
1909 asyclt(i,j)=(g1+g2)/(w1+w2)
1910
1911 endif
1912
1913 enddo
1914 enddo
1915
1916 do j=1,n
1917 do i=1,m
1918 ssacl(i,j,k)=ssaclt(i,j)
1919 enddo
1920 enddo
1921 do j=1,n
1922 do i=1,m
1923 asycl(i,j,k)=asyclt(i,j)
1924 enddo
1925 enddo
1926
1927 enddo
1928
1929 !-----integration over the k-distribution function
1930
1931 do 200 ik=1,nk
1932
1933 do 300 k= 1, np
1934
1935 do j= 1, n
1936 do i= 1, m
1937
1938 tauwv=xk(ik)*wh(i,j,k)
1939
1940 !-----compute clear-sky optical thickness, single scattering albedo,
1941 ! and asymmetry factor.
1942
1943 tausto=tauwv+taual(i,j,k,iv)+1.0e-8
1944 ssatau=ssaal(i,j,k,iv)*taual(i,j,k,iv)
1945 asysto=asyal(i,j,k,iv)*ssaal(i,j,k,iv)*taual(i,j,k,iv)
1946
1947 !-----compute reflectance and transmittance for cloudless layers
1948
1949 tauto=tausto
1950 ssato=ssatau/tauto+1.0e-8
1951
1952 if (ssato .gt. 0.001) then
1953
1954 ssato=min(ssato,0.999999)
1955 asyto=asysto/(ssato*tauto)
1956
1957 !- for direct incident radiation
1958
1959 call deledd (tauto,ssato,asyto,csm(i,j), &
1960 rr1t(i,j),tt1t(i,j),td1t(i,j))
1961
1962 !- for diffuse incident radiation
1963
1964 call sagpol (tauto,ssato,asyto,rs1t(i,j),ts1t(i,j))
1965
1966 else
1967
1968 td1t(i,j)=exp(-tauto*csm(i,j))
1969 ts1t(i,j)=exp(-1.66*tauto)
1970 tt1t(i,j)=0.0
1971 rr1t(i,j)=0.0
1972 rs1t(i,j)=0.0
1973
1974 endif
1975
1976 !-----compute reflectance and transmittance for cloud layers
1977
1978 if (tauclb(i,j,k).lt.0.01 .or. fcld(i,j,k).lt.0.01) then
1979
1980 rr2t(i,j)=rr1t(i,j)
1981 tt2t(i,j)=tt1t(i,j)
1982 td2t(i,j)=td1t(i,j)
1983 rs2t(i,j)=rs1t(i,j)
1984 ts2t(i,j)=ts1t(i,j)
1985
1986 else
1987
1988 !- for direct incident radiation
1989
1990 tauto=tausto+tauclb(i,j,k)
1991 ssato=(ssatau+ssacl(i,j,k)*tauclb(i,j,k))/tauto+1.0e-8
1992 ssato=min(ssato,0.999999)
1993 asyto=(asysto+asycl(i,j,k)*ssacl(i,j,k)*tauclb(i,j,k))/ &
1994 (ssato*tauto)
1995
1996 call deledd (tauto,ssato,asyto,csm(i,j), &
1997 rr2t(i,j),tt2t(i,j),td2t(i,j))
1998
1999 !- for diffuse incident radiation
2000
2001 tauto=tausto+tauclf(i,j,k)
2002 ssato=(ssatau+ssacl(i,j,k)*tauclf(i,j,k))/tauto+1.0e-8
2003 ssato=min(ssato,0.999999)
2004 asyto=(asysto+asycl(i,j,k)*ssacl(i,j,k)*tauclf(i,j,k))/ &
2005 (ssato*tauto)
2006
2007 call sagpol (tauto,ssato,asyto,rs2t(i,j),ts2t(i,j))
2008
2009 endif
2010
2011 enddo
2012 enddo
2013
2014 do j=1,n
2015 do i=1,m
2016 rr(i,j,k,1)=rr1t(i,j)
2017 enddo
2018 enddo
2019 do j=1,n
2020 do i=1,m
2021 tt(i,j,k,1)=tt1t(i,j)
2022 enddo
2023 enddo
2024 do j=1,n
2025 do i=1,m
2026 td(i,j,k,1)=td1t(i,j)
2027 enddo
2028 enddo
2029 do j=1,n
2030 do i=1,m
2031 rs(i,j,k,1)=rs1t(i,j)
2032 enddo
2033 enddo
2034 do j=1,n
2035 do i=1,m
2036 ts(i,j,k,1)=ts1t(i,j)
2037 enddo
2038 enddo
2039
2040 do j=1,n
2041 do i=1,m
2042 rr(i,j,k,2)=rr2t(i,j)
2043 enddo
2044 enddo
2045 do j=1,n
2046 do i=1,m
2047 tt(i,j,k,2)=tt2t(i,j)
2048 enddo
2049 enddo
2050 do j=1,n
2051 do i=1,m
2052 td(i,j,k,2)=td2t(i,j)
2053 enddo
2054 enddo
2055 do j=1,n
2056 do i=1,m
2057 rs(i,j,k,2)=rs2t(i,j)
2058 enddo
2059 enddo
2060 do j=1,n
2061 do i=1,m
2062 ts(i,j,k,2)=ts2t(i,j)
2063 enddo
2064 enddo
2065
2066 300 continue
2067
2068 !-----flux calculations
2069
2070 call cldflx (m,n,np,ict,icb,overcast,cc,rr,tt,td,rs,ts, &
2071 fclr,fall,fallu,falld,fsdir,fsdif)
2072
2073 do k= 1, np+1
2074 do j= 1, n
2075 do i= 1, m
2076 flx(i,j,k) = flx(i,j,k)+fall(i,j,k)*hk(ib,ik)
2077 flxu(i,j,k) = flxu(i,j,k)+fallu(i,j,k)*hk(ib,ik)
2078 flxd(i,j,k) = flxd(i,j,k)+falld(i,j,k)*hk(ib,ik)
2079 enddo
2080 enddo
2081 do j= 1, n
2082 do i= 1, m
2083 flc(i,j,k) = flc(i,j,k)+fclr(i,j,k)*hk(ib,ik)
2084 enddo
2085 enddo
2086 enddo
2087
2088 !-----compute downward surface fluxes in the ir region
2089
2090 do j= 1, n
2091 do i= 1, m
2092 fdirir(i,j) = fdirir(i,j)+fsdir(i,j)*hk(ib,ik)
2093 fdifir(i,j) = fdifir(i,j)+fsdif(i,j)*hk(ib,ik)
2094 enddo
2095 enddo
2096
2097 200 continue
2098 100 continue
2099
2100 end subroutine solir
2101
2102 !********************************************************************
2103
2104 subroutine cldscale (m,n,ndim,np,cosz,fcld,taucld,ict,icb, &
2105 cc,tauclb,tauclf)
2106
2107 !********************************************************************
2108 !
2109 ! This subroutine computes the high, middle, and
2110 ! low cloud amounts and scales the cloud optical thickness.
2111 !
2112 ! To simplify calculations in a cloudy atmosphere, clouds are
2113 ! grouped into high, middle and low clouds separated by the levels
2114 ! ict and icb (level 1 is the top of the model atmosphere).
2115 !
2116 ! Within each of the three groups, clouds are assumed maximally
2117 ! overlapped, and the cloud cover (cc) of a group is the maximum
2118 ! cloud cover of all the layers in the group. The optical thickness
2119 ! (taucld) of a given layer is then scaled to new values (tauclb and
2120 ! tauclf) so that the layer reflectance corresponding to the cloud
2121 ! cover cc is the same as the original reflectance with optical
2122 ! thickness taucld and cloud cover fcld.
2123 !
2124 !---input parameters
2125 !
2126 ! number of grid intervals in zonal direction (m)
2127 ! number of grid intervals in meridional direction (n)
2128 ! maximum number of grid intervals in meridional direction (ndim)
2129 ! number of atmospheric layers (np)
2130 ! cosine of the solar zenith angle (cosz)
2131 ! fractional cloud cover (fcld)
2132 ! cloud optical thickness (taucld)
2133 ! index separating high and middle clouds (ict)
2134 ! index separating middle and low clouds (icb)
2135 !
2136 !---output parameters
2137 !
2138 ! fractional cover of high, middle, and low clouds (cc)
2139 ! scaled cloud optical thickness for beam radiation (tauclb)
2140 ! scaled cloud optical thickness for diffuse radiation (tauclf)
2141 !
2142 !********************************************************************
2143 implicit none
2144 !********************************************************************
2145
2146 !-----input parameters
2147
2148 integer m,n,ndim,np
2149 integer ict(m,ndim),icb(m,ndim)
2150 real cosz(m,ndim),fcld(m,ndim,np),taucld(m,ndim,np,2)
2151
2152 !-----output parameters
2153
2154 real cc(m,n,3),tauclb(m,n,np),tauclf(m,n,np)
2155
2156 !-----temporary variables
2157
2158 integer i,j,k,im,it,ia,kk
2159 real fm,ft,fa,xai,taux
2160
2161 !-----pre-computed table
2162
2163 integer nm,nt,na
2164 parameter (nm=11,nt=9,na=11)
2165 real dm,dt,da,t1,caib(nm,nt,na),caif(nt,na)
2166 parameter (dm=0.1,dt=0.30103,da=0.1,t1=-0.9031)
2167
2168 !-----include the pre-computed table of mcai for scaling the cloud optical
2169 ! thickness under the assumption that clouds are maximally overlapped
2170 !
2171 ! caib is for scaling the cloud optical thickness for direct radiation
2172 ! caif is for scaling the cloud optical thickness for diffuse radiation
2173
2174
2175 data ((caib(1,i,j),j=1,11),i=1,9)/ &
2176 .000,0.068,0.140,0.216,0.298,0.385,0.481,0.586,0.705,0.840,1.000, &
2177 .000,0.052,0.106,0.166,0.230,0.302,0.383,0.478,0.595,0.752,1.000, &
2178 .000,0.038,0.078,0.120,0.166,0.218,0.276,0.346,0.438,0.582,1.000, &
2179 .000,0.030,0.060,0.092,0.126,0.164,0.206,0.255,0.322,0.442,1.000, &
2180 .000,0.025,0.051,0.078,0.106,0.136,0.170,0.209,0.266,0.462,1.000, &
2181 .000,0.023,0.046,0.070,0.095,0.122,0.150,0.187,0.278,0.577,1.000, &
2182 .000,0.022,0.043,0.066,0.089,0.114,0.141,0.187,0.354,0.603,1.000, &
2183 .000,0.021,0.042,0.063,0.086,0.108,0.135,0.214,0.349,0.565,1.000, &
2184 .000,0.021,0.041,0.062,0.083,0.105,0.134,0.202,0.302,0.479,1.000/
2185 data ((caib(2,i,j),j=1,11),i=1,9)/ &
2186 .000,0.088,0.179,0.272,0.367,0.465,0.566,0.669,0.776,0.886,1.000, &
2187 .000,0.079,0.161,0.247,0.337,0.431,0.531,0.637,0.749,0.870,1.000, &
2188 .000,0.065,0.134,0.207,0.286,0.372,0.466,0.572,0.692,0.831,1.000, &
2189 .000,0.049,0.102,0.158,0.221,0.290,0.370,0.465,0.583,0.745,1.000, &
2190 .000,0.037,0.076,0.118,0.165,0.217,0.278,0.354,0.459,0.638,1.000, &
2191 .000,0.030,0.061,0.094,0.130,0.171,0.221,0.286,0.398,0.631,1.000, &
2192 .000,0.026,0.052,0.081,0.111,0.146,0.189,0.259,0.407,0.643,1.000, &
2193 .000,0.023,0.047,0.072,0.098,0.129,0.170,0.250,0.387,0.598,1.000, &
2194 .000,0.022,0.044,0.066,0.090,0.118,0.156,0.224,0.328,0.508,1.000/
2195 data ((caib(3,i,j),j=1,11),i=1,9)/ &
2196 .000,0.094,0.189,0.285,0.383,0.482,0.582,0.685,0.788,0.894,1.000, &
2197 .000,0.088,0.178,0.271,0.366,0.465,0.565,0.669,0.776,0.886,1.000, &
2198 .000,0.079,0.161,0.247,0.337,0.431,0.531,0.637,0.750,0.870,1.000, &
2199 .000,0.066,0.134,0.209,0.289,0.375,0.470,0.577,0.697,0.835,1.000, &
2200 .000,0.050,0.104,0.163,0.227,0.300,0.383,0.483,0.606,0.770,1.000, &
2201 .000,0.038,0.080,0.125,0.175,0.233,0.302,0.391,0.518,0.710,1.000, &
2202 .000,0.031,0.064,0.100,0.141,0.188,0.249,0.336,0.476,0.689,1.000, &
2203 .000,0.026,0.054,0.084,0.118,0.158,0.213,0.298,0.433,0.638,1.000, &
2204 .000,0.023,0.048,0.074,0.102,0.136,0.182,0.254,0.360,0.542,1.000/
2205 data ((caib(4,i,j),j=1,11),i=1,9)/ &
2206 .000,0.096,0.193,0.290,0.389,0.488,0.589,0.690,0.792,0.896,1.000, &
2207 .000,0.092,0.186,0.281,0.378,0.477,0.578,0.680,0.785,0.891,1.000, &
2208 .000,0.086,0.174,0.264,0.358,0.455,0.556,0.660,0.769,0.882,1.000, &
2209 .000,0.074,0.153,0.235,0.323,0.416,0.514,0.622,0.737,0.862,1.000, &
2210 .000,0.061,0.126,0.195,0.271,0.355,0.449,0.555,0.678,0.823,1.000, &
2211 .000,0.047,0.098,0.153,0.215,0.286,0.370,0.471,0.600,0.770,1.000, &
2212 .000,0.037,0.077,0.120,0.170,0.230,0.303,0.401,0.537,0.729,1.000, &
2213 .000,0.030,0.062,0.098,0.138,0.187,0.252,0.343,0.476,0.673,1.000, &
2214 .000,0.026,0.053,0.082,0.114,0.154,0.207,0.282,0.391,0.574,1.000/
2215 data ((caib(5,i,j),j=1,11),i=1,9)/ &
2216 .000,0.097,0.194,0.293,0.392,0.492,0.592,0.693,0.794,0.897,1.000, &
2217 .000,0.094,0.190,0.286,0.384,0.483,0.584,0.686,0.789,0.894,1.000, &
2218 .000,0.090,0.181,0.274,0.370,0.468,0.569,0.672,0.778,0.887,1.000, &
2219 .000,0.081,0.165,0.252,0.343,0.439,0.539,0.645,0.757,0.874,1.000, &
2220 .000,0.069,0.142,0.218,0.302,0.392,0.490,0.598,0.717,0.850,1.000, &
2221 .000,0.054,0.114,0.178,0.250,0.330,0.422,0.529,0.656,0.810,1.000, &
2222 .000,0.042,0.090,0.141,0.200,0.269,0.351,0.455,0.589,0.764,1.000, &
2223 .000,0.034,0.070,0.112,0.159,0.217,0.289,0.384,0.515,0.703,1.000, &
2224 .000,0.028,0.058,0.090,0.128,0.174,0.231,0.309,0.420,0.602,1.000/
2225 data ((caib(6,i,j),j=1,11),i=1,9)/ &
2226 .000,0.098,0.196,0.295,0.394,0.494,0.594,0.695,0.796,0.898,1.000, &
2227 .000,0.096,0.193,0.290,0.389,0.488,0.588,0.690,0.792,0.895,1.000, &
2228 .000,0.092,0.186,0.281,0.378,0.477,0.577,0.680,0.784,0.891,1.000, &
2229 .000,0.086,0.174,0.264,0.358,0.455,0.556,0.661,0.769,0.882,1.000, &
2230 .000,0.075,0.154,0.237,0.325,0.419,0.518,0.626,0.741,0.865,1.000, &
2231 .000,0.062,0.129,0.201,0.279,0.366,0.462,0.571,0.694,0.836,1.000, &
2232 .000,0.049,0.102,0.162,0.229,0.305,0.394,0.501,0.631,0.793,1.000, &
2233 .000,0.038,0.080,0.127,0.182,0.245,0.323,0.422,0.550,0.730,1.000, &
2234 .000,0.030,0.064,0.100,0.142,0.192,0.254,0.334,0.448,0.627,1.000/
2235 data ((caib(7,i,j),j=1,11),i=1,9)/ &
2236 .000,0.098,0.198,0.296,0.396,0.496,0.596,0.696,0.797,0.898,1.000, &
2237 .000,0.097,0.194,0.293,0.392,0.491,0.591,0.693,0.794,0.897,1.000, &
2238 .000,0.094,0.190,0.286,0.384,0.483,0.583,0.686,0.789,0.894,1.000, &
2239 .000,0.089,0.180,0.274,0.369,0.467,0.568,0.672,0.778,0.887,1.000, &
2240 .000,0.081,0.165,0.252,0.344,0.440,0.541,0.646,0.758,0.875,1.000, &
2241 .000,0.069,0.142,0.221,0.306,0.397,0.496,0.604,0.722,0.854,1.000, &
2242 .000,0.056,0.116,0.182,0.256,0.338,0.432,0.540,0.666,0.816,1.000, &
2243 .000,0.043,0.090,0.143,0.203,0.273,0.355,0.455,0.583,0.754,1.000, &
2244 .000,0.034,0.070,0.111,0.157,0.210,0.276,0.359,0.474,0.650,1.000/
2245 data ((caib(8,i,j),j=1,11),i=1,9)/ &
2246 .000,0.099,0.198,0.298,0.398,0.497,0.598,0.698,0.798,0.899,1.000, &
2247 .000,0.098,0.196,0.295,0.394,0.494,0.594,0.695,0.796,0.898,1.000, &
2248 .000,0.096,0.193,0.290,0.390,0.489,0.589,0.690,0.793,0.896,1.000, &
2249 .000,0.093,0.186,0.282,0.379,0.478,0.578,0.681,0.786,0.892,1.000, &
2250 .000,0.086,0.175,0.266,0.361,0.458,0.558,0.663,0.771,0.883,1.000, &
2251 .000,0.076,0.156,0.240,0.330,0.423,0.523,0.630,0.744,0.867,1.000, &
2252 .000,0.063,0.130,0.203,0.282,0.369,0.465,0.572,0.694,0.834,1.000, &
2253 .000,0.049,0.102,0.161,0.226,0.299,0.385,0.486,0.611,0.774,1.000, &
2254 .000,0.038,0.078,0.122,0.172,0.229,0.297,0.382,0.498,0.672,1.000/
2255 data ((caib(9,i,j),j=1,11),i=1,9)/ &
2256 .000,0.099,0.199,0.298,0.398,0.498,0.598,0.699,0.799,0.899,1.000, &
2257 .000,0.099,0.198,0.298,0.398,0.497,0.598,0.698,0.798,0.899,1.000, &
2258 .000,0.098,0.196,0.295,0.394,0.494,0.594,0.695,0.796,0.898,1.000, &
2259 .000,0.096,0.193,0.290,0.389,0.488,0.588,0.690,0.792,0.895,1.000, &
2260 .000,0.092,0.185,0.280,0.376,0.474,0.575,0.678,0.782,0.890,1.000, &
2261 .000,0.084,0.170,0.259,0.351,0.447,0.547,0.652,0.762,0.878,1.000, &
2262 .000,0.071,0.146,0.224,0.308,0.398,0.494,0.601,0.718,0.850,1.000, &
2263 .000,0.056,0.114,0.178,0.248,0.325,0.412,0.514,0.638,0.793,1.000, &
2264 .000,0.042,0.086,0.134,0.186,0.246,0.318,0.405,0.521,0.691,1.000/
2265 data ((caib(10,i,j),j=1,11),i=1,9)/ &
2266 .000,0.100,0.200,0.300,0.400,0.500,0.600,0.700,0.800,0.900,1.000, &
2267 .000,0.100,0.200,0.300,0.400,0.500,0.600,0.700,0.800,0.900,1.000, &
2268 .000,0.100,0.200,0.300,0.400,0.500,0.600,0.700,0.800,0.900,1.000, &
2269 .000,0.100,0.199,0.298,0.398,0.498,0.598,0.698,0.798,0.899,1.000, &
2270 .000,0.098,0.196,0.294,0.392,0.491,0.590,0.691,0.793,0.896,1.000, &
2271 .000,0.092,0.185,0.278,0.374,0.470,0.570,0.671,0.777,0.886,1.000, &
2272 .000,0.081,0.162,0.246,0.333,0.424,0.521,0.625,0.738,0.862,1.000, &
2273 .000,0.063,0.128,0.196,0.270,0.349,0.438,0.540,0.661,0.809,1.000, &
2274 .000,0.046,0.094,0.146,0.202,0.264,0.337,0.426,0.542,0.710,1.000/
2275 data ((caib(11,i,j),j=1,11),i=1,9)/ &
2276 .000,0.101,0.202,0.302,0.402,0.502,0.602,0.702,0.802,0.901,1.000, &
2277 .000,0.102,0.202,0.303,0.404,0.504,0.604,0.703,0.802,0.902,1.000, &
2278 .000,0.102,0.205,0.306,0.406,0.506,0.606,0.706,0.804,0.902,1.000, &
2279 .000,0.104,0.207,0.309,0.410,0.510,0.609,0.707,0.805,0.902,1.000, &
2280 .000,0.106,0.208,0.309,0.409,0.508,0.606,0.705,0.803,0.902,1.000, &
2281 .000,0.102,0.202,0.298,0.395,0.493,0.590,0.690,0.790,0.894,1.000, &
2282 .000,0.091,0.179,0.267,0.357,0.449,0.545,0.647,0.755,0.872,1.000, &
2283 .000,0.073,0.142,0.214,0.290,0.372,0.462,0.563,0.681,0.822,1.000, &
2284 .000,0.053,0.104,0.158,0.217,0.281,0.356,0.446,0.562,0.726,1.000/
2285 data ((caif(i,j),j=1,11),i=1,9)/ &
2286 .000,0.099,0.198,0.297,0.397,0.496,0.597,0.697,0.798,0.899,1.000, &
2287 .000,0.098,0.196,0.294,0.394,0.494,0.594,0.694,0.796,0.898,1.000, &
2288 .000,0.096,0.192,0.290,0.388,0.487,0.587,0.689,0.792,0.895,1.000, &
2289 .000,0.092,0.185,0.280,0.376,0.476,0.576,0.678,0.783,0.890,1.000, &
2290 .000,0.085,0.173,0.263,0.357,0.454,0.555,0.659,0.768,0.881,1.000, &
2291 .000,0.076,0.154,0.237,0.324,0.418,0.517,0.624,0.738,0.864,1.000, &
2292 .000,0.063,0.131,0.203,0.281,0.366,0.461,0.567,0.688,0.830,1.000, &
2293 .000,0.052,0.107,0.166,0.232,0.305,0.389,0.488,0.610,0.770,1.000, &
2294 .000,0.043,0.088,0.136,0.189,0.248,0.317,0.400,0.510,0.675,1.000/
2295
2296 !-----clouds within each of the high, middle, and low clouds are assumed
2297 ! to be maximally overlapped, and the cloud cover (cc) for a group
2298 ! (high, middle, or low) is the maximum cloud cover of all the layers
2299 ! within a group
2300
2301 do j=1,n
2302 do i=1,m
2303 cc(i,j,1)=0.0
2304 cc(i,j,2)=0.0
2305 cc(i,j,3)=0.0
2306 enddo
2307 enddo
2308 do j=1,n
2309 do i=1,m
2310 do k=1,ict(i,j)-1
2311 cc(i,j,1)=max(cc(i,j,1),fcld(i,j,k))
2312 enddo
2313 enddo
2314 enddo
2315
2316 do j=1,n
2317 do i=1,m
2318 do k=ict(i,j),icb(i,j)-1
2319 cc(i,j,2)=max(cc(i,j,2),fcld(i,j,k))
2320 enddo
2321 enddo
2322 enddo
2323
2324 do j=1,n
2325 do i=1,m
2326 do k=icb(i,j),np
2327 cc(i,j,3)=max(cc(i,j,3),fcld(i,j,k))
2328 enddo
2329 enddo
2330 enddo
2331
2332 !-----scale the cloud optical thickness.
2333 ! taucld(i,j,k,1) is the optical thickness for ice particles, and
2334 ! taucld(i,j,k,2) is the optical thickness for liquid particles.
2335
2336 do j=1,n
2337 do i=1,m
2338
2339 do k=1,np
2340
2341 if(k.lt.ict(i,j)) then
2342 kk=1
2343 elseif(k.ge.ict(i,j) .and. k.lt.icb(i,j)) then
2344 kk=2
2345 else
2346 kk=3
2347 endif
2348
2349 tauclb(i,j,k) = 0.0
2350 tauclf(i,j,k) = 0.0
2351
2352 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
2353 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
2354
2355 !-----normalize cloud cover
2356
2357 fa=fcld(i,j,k)/cc(i,j,kk)
2358
2359 !-----table look-up
2360
2361 taux=min(taux,32.)
2362
2363 fm=cosz(i,j)/dm
2364 ft=(log10(taux)-t1)/dt
2365 fa=fa/da
2366
2367 im=int(fm+1.5)
2368 it=int(ft+1.5)
2369 ia=int(fa+1.5)
2370
2371 im=max(im,2)
2372 it=max(it,2)
2373 ia=max(ia,2)
2374
2375 im=min(im,nm-1)
2376 it=min(it,nt-1)
2377 ia=min(ia,na-1)
2378
2379 fm=fm-float(im-1)
2380 ft=ft-float(it-1)
2381 fa=fa-float(ia-1)
2382
2383 !-----scale cloud optical thickness for beam radiation.
2384 ! the scaling factor, xai, is a function of the solar zenith
2385 ! angle, optical thickness, and cloud cover.
2386
2387 xai= (-caib(im-1,it,ia)*(1.-fm)+ &
2388 caib(im+1,it,ia)*(1.+fm))*fm*.5+caib(im,it,ia)*(1.-fm*fm)
2389
2390 xai=xai+(-caib(im,it-1,ia)*(1.-ft)+ &
2391 caib(im,it+1,ia)*(1.+ft))*ft*.5+caib(im,it,ia)*(1.-ft*ft)
2392
2393 xai=xai+(-caib(im,it,ia-1)*(1.-fa)+ &
2394 caib(im,it,ia+1)*(1.+fa))*fa*.5+caib(im,it,ia)*(1.-fa*fa)
2395
2396 xai= xai-2.*caib(im,it,ia)
2397 xai=max(xai,0.0)
2398
2399 tauclb(i,j,k) = taux*xai
2400
2401 !-----scale cloud optical thickness for diffuse radiation.
2402 ! the scaling factor, xai, is a function of the cloud optical
2403 ! thickness and cover but not the solar zenith angle.
2404
2405 xai= (-caif(it-1,ia)*(1.-ft)+ &
2406 caif(it+1,ia)*(1.+ft))*ft*.5+caif(it,ia)*(1.-ft*ft)
2407
2408 xai=xai+(-caif(it,ia-1)*(1.-fa)+ &
2409 caif(it,ia+1)*(1.+fa))*fa*.5+caif(it,ia)*(1.-fa*fa)
2410
2411 xai= xai-caif(it,ia)
2412 xai=max(xai,0.0)
2413
2414 tauclf(i,j,k) = taux*xai
2415
2416 endif
2417
2418 enddo
2419 enddo
2420 enddo
2421
2422 end subroutine cldscale
2423
2424 !*********************************************************************
2425
2426 subroutine deledd(tau,ssc,g0,csm,rr,tt,td)
2427
2428 !*********************************************************************
2429 !
2430 !-----uses the delta-eddington approximation to compute the
2431 ! bulk scattering properties of a single layer
2432 ! coded following King and Harshvardhan (JAS, 1986)
2433 !
2434 ! inputs:
2435 !
2436 ! tau: the effective optical thickness
2437 ! ssc: the effective single scattering albedo
2438 ! g0: the effective asymmetry factor
2439 ! csm: the effective secant of the zenith angle
2440 !
2441 ! outputs:
2442 !
2443 ! rr: the layer reflection of the direct beam
2444 ! tt: the layer diffuse transmission of the direct beam
2445 ! td: the layer direct transmission of the direct beam
2446 !
2447 !*********************************************************************
2448 implicit none
2449 !*********************************************************************
2450
2451 real zero,one,two,three,four,fourth,seven,thresh
2452 parameter (one =1., three=3.)
2453 parameter (two =2., seven=7.)
2454 parameter (four=4., fourth=.25)
2455 parameter (zero=0., thresh=1.e-8)
2456
2457 !-----input parameters
2458 real tau,ssc,g0,csm
2459
2460 !-----output parameters
2461 real rr,tt,td
2462
2463 !-----temporary parameters
2464
2465 real zth,ff,xx,taup,sscp,gp,gm1,gm2,gm3,akk,alf1,alf2, &
2466 all,bll,st7,st8,cll,dll,fll,ell,st1,st2,st3,st4
2467
2468 !---------------------------------------------------------------------
2469
2470 zth = one / csm
2471
2472 ! delta-eddington scaling of single scattering albedo,
2473 ! optical thickness, and asymmetry factor,
2474 ! K & H eqs(27-29)
2475
2476 ff = g0*g0
2477 xx = one-ff*ssc
2478 taup= tau*xx
2479 sscp= ssc*(one-ff)/xx
2480 gp = g0/(one+g0)
2481
2482 ! gamma1, gamma2, and gamma3. see table 2 and eq(26) K & H
2483 ! ssc and gp are the d-s single scattering
2484 ! albedo and asymmetry factor.
2485
2486 xx = three*gp
2487 gm1 = (seven - sscp*(four+xx))*fourth
2488 gm2 = -(one - sscp*(four-xx))*fourth
2489
2490 ! akk is k as defined in eq(25) of K & H
2491
2492 akk = sqrt((gm1+gm2)*(gm1-gm2))
2493
2494 xx = akk * zth
2495 st7 = one - xx
2496 st8 = one + xx
2497 st3 = st7 * st8
2498
2499 if (abs(st3) .lt. thresh) then
2500 zth = zth + 0.001
2501 xx = akk * zth
2502 st7 = one - xx
2503 st8 = one + xx
2504 st3 = st7 * st8
2505 endif
2506
2507 ! extinction of the direct beam transmission
2508
2509 td = exp(-taup/zth)
2510
2511 ! alf1 and alf2 are alpha1 and alpha2 from eqs (23) & (24) of K & H
2512
2513 gm3 = (two - zth*three*gp)*fourth
2514 xx = gm1 - gm2
2515 alf1 = gm1 - gm3 * xx
2516 alf2 = gm2 + gm3 * xx
2517
2518 ! all is last term in eq(21) of K & H
2519 ! bll is last term in eq(22) of K & H
2520
2521 xx = akk * two
2522 all = (gm3 - alf2 * zth )*xx*td
2523 bll = (one - gm3 + alf1*zth)*xx
2524
2525 xx = akk * gm3
2526 cll = (alf2 + xx) * st7
2527 dll = (alf2 - xx) * st8
2528
2529 xx = akk * (one-gm3)
2530 fll = (alf1 + xx) * st8
2531 ell = (alf1 - xx) * st7
2532
2533 st2 = exp(-akk*taup)
2534 st4 = st2 * st2
2535
2536 st1 = sscp / ((akk+gm1 + (akk-gm1)*st4) * st3)
2537
2538 ! rr is r-hat of eq(21) of K & H
2539 ! tt is diffuse part of t-hat of eq(22) of K & H
2540
2541 rr = ( cll-dll*st4 -all*st2)*st1
2542 tt = - ((fll-ell*st4)*td-bll*st2)*st1
2543
2544 rr = max(rr,zero)
2545 tt = max(tt,zero)
2546
2547 end subroutine deledd
2548
2549 !*********************************************************************
2550
2551 subroutine sagpol(tau,ssc,g0,rll,tll)
2552
2553 !*********************************************************************
2554 !-----transmittance (tll) and reflectance (rll) of diffuse radiation
2555 ! follows Sagan and Pollock (JGR, 1967).
2556 ! also, eq.(31) of Lacis and Hansen (JAS, 1974).
2557 !
2558 !-----input parameters:
2559 !
2560 ! tau: the effective optical thickness
2561 ! ssc: the effective single scattering albedo
2562 ! g0: the effective asymmetry factor
2563 !
2564 !-----output parameters:
2565 !
2566 ! rll: the layer reflection of diffuse radiation
2567 ! tll: the layer transmission of diffuse radiation
2568 !
2569 !*********************************************************************
2570 implicit none
2571 !*********************************************************************
2572
2573 real one,three,four
2574 parameter (one=1., three=3., four=4.)
2575
2576 !-----output parameters:
2577
2578 real tau,ssc,g0
2579
2580 !-----output parameters:
2581
2582 real rll,tll
2583
2584 !-----temporary arrays
2585
2586 real xx,uuu,ttt,emt,up1,um1,st1
2587
2588 xx = one-ssc*g0
2589 uuu = sqrt( xx/(one-ssc))
2590 ttt = sqrt( xx*(one-ssc)*three )*tau
2591 emt = exp(-ttt)
2592 up1 = uuu + one
2593 um1 = uuu - one
2594 xx = um1*emt
2595 st1 = one / ((up1+xx) * (up1-xx))
2596 rll = up1*um1*(one-emt*emt)*st1
2597 tll = uuu*four*emt *st1
2598
2599 end subroutine sagpol
2600
2601 !*******************************************************************
2602
2603 subroutine cldflx (m,n,np,ict,icb,overcast,cc,rr,tt,td,rs,ts,&
2604 fclr,fall,fallu,falld,fsdir,fsdif)
2605
2606 !*******************************************************************
2607 ! compute upward and downward fluxes using a two-stream adding method
2608 ! following equations (3)-(5) of Chou (1992, JAS).
2609 !
2610 ! clouds are grouped into high, middle, and low clouds which are
2611 ! assumed randomly overlapped. It involves eight sets of calculations.
2612 ! In each set of calculations, each atmospheric layer is homogeneous,
2613 ! either totally filled with clouds or without clouds.
2614
2615 ! input parameters:
2616 !
2617 ! m: number of soundings in zonal direction
2618 ! n: number of soundings in meridional direction
2619 ! np: number of atmospheric layers
2620 ! ict: the level separating high and middle clouds
2621 ! icb: the level separating middle and low clouds
2622 ! cc: effective cloud covers for high, middle and low clouds
2623 ! tt: diffuse transmission of a layer illuminated by beam radiation
2624 ! td: direct beam tranmssion
2625 ! ts: transmission of a layer illuminated by diffuse radiation
2626 ! rr: reflection of a layer illuminated by beam radiation
2627 ! rs: reflection of a layer illuminated by diffuse radiation
2628 !
2629 ! output parameters:
2630 !
2631 ! fclr: clear-sky flux (downward minus upward)
2632 ! fall: all-sky flux (downward minus upward)
2633 ! fsdir: surface direct downward flux
2634 ! fsdif: surface diffuse downward flux
2635 !
2636 !*********************************************************************c
2637 implicit none
2638 !*********************************************************************c
2639
2640 !-----input parameters
2641
2642 integer m,n,np
2643 integer ict(m,n),icb(m,n)
2644
2645 real rr(m,n,np+1,2),tt(m,n,np+1,2),td(m,n,np+1,2)
2646 real rs(m,n,np+1,2),ts(m,n,np+1,2)
2647 real cc(m,n,3)
2648 logical overcast
2649
2650 !-----temporary array
2651
2652 integer i,j,k,ih,im,is,itm
2653 real rra(m,n,np+1,2,2),tta(m,n,np+1,2,2),tda(m,n,np+1,2,2)
2654 real rsa(m,n,np+1,2,2),rxa(m,n,np+1,2,2)
2655 real ch(m,n),cm(m,n),ct(m,n),flxdn(m,n,np+1)
2656 real flxdnu(m,n,np+1),flxdnd(m,n,np+1)
2657 real fdndir(m,n),fdndif(m,n),fupdif
2658 real denm,xx
2659
2660 !-----output parameters
2661
2662 real fclr(m,n,np+1),fall(m,n,np+1)
2663 real fallu(m,n,np+1),falld(m,n,np+1)
2664 real fsdir(m,n),fsdif(m,n)
2665
2666 !-----initialize all-sky flux (fall) and surface downward fluxes
2667
2668 do k=1,np+1
2669 do j=1,n
2670 do i=1,m
2671 fclr(i,j,k)=0.0
2672 fall(i,j,k)=0.0
2673 fallu(i,j,k)=0.0
2674 falld(i,j,k)=0.0
2675 enddo
2676 enddo
2677 enddo
2678
2679 do j=1,n
2680 do i=1,m
2681 fsdir(i,j)=0.0
2682 fsdif(i,j)=0.0
2683 enddo
2684 enddo
2685
2686 !-----compute transmittances and reflectances for a composite of
2687 ! layers. layers are added one at a time, going down from the top.
2688 ! tda is the composite transmittance illuminated by beam radiation
2689 ! tta is the composite diffuse transmittance illuminated by
2690 ! beam radiation
2691 ! rsa is the composite reflectance illuminated from below
2692 ! by diffuse radiation
2693 ! tta and rsa are computed from eqs. (4b) and (3b) of Chou
2694
2695 itm=1
2696
2697 !-----if overcas.=.true., set itm=2, and only one set of fluxes is computed
2698
2699 if (overcast) itm=2
2700
2701 !-----for high clouds. indices 1 and 2 denote clear and cloudy
2702 ! situations, respectively.
2703
2704 do 10 ih=itm,2
2705
2706 do j= 1, n
2707 do i= 1, m
2708 tda(i,j,1,ih,1)=td(i,j,1,ih)
2709 tta(i,j,1,ih,1)=tt(i,j,1,ih)
2710 rsa(i,j,1,ih,1)=rs(i,j,1,ih)
2711 tda(i,j,1,ih,2)=td(i,j,1,ih)
2712 tta(i,j,1,ih,2)=tt(i,j,1,ih)
2713 rsa(i,j,1,ih,2)=rs(i,j,1,ih)
2714 enddo
2715 enddo
2716
2717 do j= 1, n
2718 do i= 1, m
2719 do k= 2, ict(i,j)-1
2720 denm = ts(i,j,k,ih)/( 1.-rsa(i,j,k-1,ih,1)*rs(i,j,k,ih))
2721 tda(i,j,k,ih,1)= tda(i,j,k-1,ih,1)*td(i,j,k,ih)
2722 tta(i,j,k,ih,1)= tda(i,j,k-1,ih,1)*tt(i,j,k,ih) &
2723 +(tda(i,j,k-1,ih,1)*rr(i,j,k,ih) &
2724 *rsa(i,j,k-1,ih,1)+tta(i,j,k-1,ih,1))*denm
2725 rsa(i,j,k,ih,1)= rs(i,j,k,ih)+ts(i,j,k,ih) &
2726 *rsa(i,j,k-1,ih,1)*denm
2727 tda(i,j,k,ih,2)= tda(i,j,k,ih,1)
2728 tta(i,j,k,ih,2)= tta(i,j,k,ih,1)
2729 rsa(i,j,k,ih,2)= rsa(i,j,k,ih,1)
2730 enddo
2731 enddo
2732 enddo
2733
2734 !-----for middle clouds
2735
2736 do 10 im=itm,2
2737
2738 do j= 1, n
2739 do i= 1, m
2740 do k= ict(i,j), icb(i,j)-1
2741 denm = ts(i,j,k,im)/( 1.-rsa(i,j,k-1,ih,im)*rs(i,j,k,im))
2742 tda(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*td(i,j,k,im)
2743 tta(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*tt(i,j,k,im) &
2744 +(tda(i,j,k-1,ih,im)*rr(i,j,k,im) &
2745 *rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2746 rsa(i,j,k,ih,im)= rs(i,j,k,im)+ts(i,j,k,im) &
2747 *rsa(i,j,k-1,ih,im)*denm
2748 enddo
2749 enddo
2750 enddo
2751
2752 10 continue
2753
2754 !-----layers are added one at a time, going up from the surface.
2755 ! rra is the composite reflectance illuminated by beam radiation
2756 ! rxa is the composite reflectance illuminated from above
2757 ! by diffuse radiation
2758 ! rra and rxa are computed from eqs. (4a) and (3a) of Chou
2759
2760 !-----for the low clouds
2761
2762 do 20 is=itm,2
2763
2764 do j= 1, n
2765 do i= 1, m
2766 rra(i,j,np+1,1,is)=rr(i,j,np+1,is)
2767 rxa(i,j,np+1,1,is)=rs(i,j,np+1,is)
2768 rra(i,j,np+1,2,is)=rr(i,j,np+1,is)
2769 rxa(i,j,np+1,2,is)=rs(i,j,np+1,is)
2770 enddo
2771 enddo
2772
2773 do j= 1, n
2774 do i= 1, m
2775 do k=np,icb(i,j),-1
2776 denm=ts(i,j,k,is)/( 1.-rs(i,j,k,is)*rxa(i,j,k+1,1,is) )
2777 rra(i,j,k,1,is)=rr(i,j,k,is)+(td(i,j,k,is) &
2778 *rra(i,j,k+1,1,is)+tt(i,j,k,is)*rxa(i,j,k+1,1,is))*denm
2779 rxa(i,j,k,1,is)= rs(i,j,k,is)+ts(i,j,k,is) &
2780 *rxa(i,j,k+1,1,is)*denm
2781 rra(i,j,k,2,is)=rra(i,j,k,1,is)
2782 rxa(i,j,k,2,is)=rxa(i,j,k,1,is)
2783 enddo
2784 enddo
2785 enddo
2786
2787 !-----for middle clouds
2788
2789 do 20 im=itm,2
2790
2791 do j= 1, n
2792 do i= 1, m
2793 do k= icb(i,j)-1,ict(i,j),-1
2794 denm=ts(i,j,k,im)/( 1.-rs(i,j,k,im)*rxa(i,j,k+1,im,is) )
2795 rra(i,j,k,im,is)= rr(i,j,k,im)+(td(i,j,k,im) &
2796 *rra(i,j,k+1,im,is)+tt(i,j,k,im)*rxa(i,j,k+1,im,is))*denm
2797 rxa(i,j,k,im,is)= rs(i,j,k,im)+ts(i,j,k,im) &
2798 *rxa(i,j,k+1,im,is)*denm
2799 enddo
2800 enddo
2801 enddo
2802
2803 20 continue
2804
2805 !-----integration over eight sky situations.
2806 ! ih, im, is denotes high, middle and low cloud groups.
2807
2808 do 100 ih=itm,2
2809
2810 !-----clear portion
2811
2812 if(ih.eq.1) then
2813 do j=1,n
2814 do i=1,m
2815 ch(i,j)=1.0-cc(i,j,1)
2816 enddo
2817 enddo
2818
2819 else
2820
2821 !-----cloudy portion
2822
2823 do j=1,n
2824 do i=1,m
2825 ch(i,j)=cc(i,j,1)
2826 enddo
2827 enddo
2828
2829 endif
2830
2831 do 100 im=itm,2
2832
2833 !-----clear portion
2834
2835 if(im.eq.1) then
2836
2837 do j=1,n
2838 do i=1,m
2839 cm(i,j)=ch(i,j)*(1.0-cc(i,j,2))
2840 enddo
2841 enddo
2842
2843 else
2844
2845 !-----cloudy portion
2846
2847 do j=1,n
2848 do i=1,m
2849 cm(i,j)=ch(i,j)*cc(i,j,2)
2850 enddo
2851 enddo
2852
2853 endif
2854
2855 do 100 is=itm,2
2856
2857 !-----clear portion
2858
2859 if(is.eq.1) then
2860
2861 do j=1,n
2862 do i=1,m
2863 ct(i,j)=cm(i,j)*(1.0-cc(i,j,3))
2864 enddo
2865 enddo
2866
2867 else
2868
2869 !-----cloudy portion
2870
2871 do j=1,n
2872 do i=1,m
2873 ct(i,j)=cm(i,j)*cc(i,j,3)
2874 enddo
2875 enddo
2876
2877 endif
2878
2879 !-----add one layer at a time, going down.
2880
2881 do j= 1, n
2882 do i= 1, m
2883 do k= icb(i,j), np
2884 denm = ts(i,j,k,is)/( 1.-rsa(i,j,k-1,ih,im)*rs(i,j,k,is) )
2885 tda(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*td(i,j,k,is)
2886 tta(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*tt(i,j,k,is) &
2887 +(tda(i,j,k-1,ih,im)*rr(i,j,k,is) &
2888 *rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2889 rsa(i,j,k,ih,im)= rs(i,j,k,is)+ts(i,j,k,is) &
2890 *rsa(i,j,k-1,ih,im)*denm
2891 enddo
2892 enddo
2893 enddo
2894
2895 !-----add one layer at a time, going up.
2896
2897 do j= 1, n
2898 do i= 1, m
2899 do k= ict(i,j)-1,1,-1
2900 denm =ts(i,j,k,ih)/(1.-rs(i,j,k,ih)*rxa(i,j,k+1,im,is))
2901 rra(i,j,k,im,is)= rr(i,j,k,ih)+(td(i,j,k,ih) &
2902 *rra(i,j,k+1,im,is)+tt(i,j,k,ih)*rxa(i,j,k+1,im,is))*denm
2903 rxa(i,j,k,im,is)= rs(i,j,k,ih)+ts(i,j,k,ih) &
2904 *rxa(i,j,k+1,im,is)*denm
2905 enddo
2906 enddo
2907 enddo
2908
2909 !-----compute fluxes following eq (5) of Chou (1992)
2910
2911 ! fdndir is the direct downward flux
2912 ! fdndif is the diffuse downward flux
2913 ! fupdif is the diffuse upward flux
2914
2915 do k=2,np+1
2916 do j=1, n
2917 do i=1, m
2918 denm= 1./(1.- rxa(i,j,k,im,is)*rsa(i,j,k-1,ih,im))
2919 fdndir(i,j)= tda(i,j,k-1,ih,im)
2920 xx = tda(i,j,k-1,ih,im)*rra(i,j,k,im,is)
2921 fdndif(i,j)= (xx*rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2922 fupdif= (xx+tta(i,j,k-1,ih,im)*rxa(i,j,k,im,is))*denm
2923 flxdn(i,j,k)=fdndir(i,j)+fdndif(i,j)-fupdif
2924 flxdnu(i,j,k)=-fupdif
2925 flxdnd(i,j,k)=fdndir(i,j)+fdndif(i,j)
2926 enddo
2927 enddo
2928 enddo
2929
2930 do j=1, n
2931 do i=1, m
2932 flxdn(i,j,1)=1.0-rra(i,j,1,im,is)
2933 flxdnu(i,j,1)=-rra(i,j,1,im,is)
2934 flxdnd(i,j,1)=1.0
2935 enddo
2936 enddo
2937
2938 !-----summation of fluxes over all (eight) sky situations.
2939
2940 do k=1,np+1
2941 do j=1,n
2942 do i=1,m
2943 if(ih.eq.1 .and. im.eq.1 .and. is.eq.1) then
2944 fclr(i,j,k)=flxdn(i,j,k)
2945 endif
2946 fall(i,j,k)=fall(i,j,k)+flxdn(i,j,k)*ct(i,j)
2947 fallu(i,j,k)=fallu(i,j,k)+flxdnu(i,j,k)*ct(i,j)
2948 falld(i,j,k)=falld(i,j,k)+flxdnd(i,j,k)*ct(i,j)
2949 enddo
2950 enddo
2951 enddo
2952
2953 do j=1,n
2954 do i=1,m
2955 fsdir(i,j)=fsdir(i,j)+fdndir(i,j)*ct(i,j)
2956 fsdif(i,j)=fsdif(i,j)+fdndif(i,j)*ct(i,j)
2957 enddo
2958 enddo
2959
2960 100 continue
2961
2962 end subroutine cldflx
2963
2964 !*****************************************************************
2965
2966 subroutine flxco2(m,n,np,swc,swh,csm,df)
2967
2968 !*****************************************************************
2969
2970 !-----compute the reduction of clear-sky downward solar flux
2971 ! due to co2 absorption.
2972
2973 implicit none
2974
2975 !-----input parameters
2976
2977 integer m,n,np
2978 real csm(m,n),swc(m,n,np+1),swh(m,n,np+1),cah(22,19)
2979
2980 !-----output (undated) parameter
2981
2982 real df(m,n,np+1)
2983
2984 !-----temporary array
2985
2986 integer i,j,k,ic,iw
2987 real xx,clog,wlog,dc,dw,x1,x2,y2
2988
2989 !********************************************************************
2990 !-----include co2 look-up table
2991
2992 data ((cah(i,j),i=1,22),j= 1, 5)/ &
2993 0.9923, 0.9922, 0.9921, 0.9920, 0.9916, 0.9910, 0.9899, 0.9882, &
2994 0.9856, 0.9818, 0.9761, 0.9678, 0.9558, 0.9395, 0.9188, 0.8945, &
2995 0.8675, 0.8376, 0.8029, 0.7621, 0.7154, 0.6647, 0.9876, 0.9876, &
2996 0.9875, 0.9873, 0.9870, 0.9864, 0.9854, 0.9837, 0.9811, 0.9773, &
2997 0.9718, 0.9636, 0.9518, 0.9358, 0.9153, 0.8913, 0.8647, 0.8350, &
2998 0.8005, 0.7599, 0.7133, 0.6627, 0.9808, 0.9807, 0.9806, 0.9805, &
2999 0.9802, 0.9796, 0.9786, 0.9769, 0.9744, 0.9707, 0.9653, 0.9573, &
3000 0.9459, 0.9302, 0.9102, 0.8866, 0.8604, 0.8311, 0.7969, 0.7565, &
3001 0.7101, 0.6596, 0.9708, 0.9708, 0.9707, 0.9705, 0.9702, 0.9697, &
3002 0.9687, 0.9671, 0.9647, 0.9612, 0.9560, 0.9483, 0.9372, 0.9221, &
3003 0.9027, 0.8798, 0.8542, 0.8253, 0.7916, 0.7515, 0.7054, 0.6551, &
3004 0.9568, 0.9568, 0.9567, 0.9565, 0.9562, 0.9557, 0.9548, 0.9533, &
3005 0.9510, 0.9477, 0.9428, 0.9355, 0.9250, 0.9106, 0.8921, 0.8700, &
3006 0.8452, 0.8171, 0.7839, 0.7443, 0.6986, 0.6486/
3007
3008 data ((cah(i,j),i=1,22),j= 6,10)/ &
3009 0.9377, 0.9377, 0.9376, 0.9375, 0.9372, 0.9367, 0.9359, 0.9345, &
3010 0.9324, 0.9294, 0.9248, 0.9181, 0.9083, 0.8948, 0.8774, 0.8565, &
3011 0.8328, 0.8055, 0.7731, 0.7342, 0.6890, 0.6395, 0.9126, 0.9126, &
3012 0.9125, 0.9124, 0.9121, 0.9117, 0.9110, 0.9098, 0.9079, 0.9052, &
3013 0.9012, 0.8951, 0.8862, 0.8739, 0.8579, 0.8385, 0.8161, 0.7900, &
3014 0.7585, 0.7205, 0.6760, 0.6270, 0.8809, 0.8809, 0.8808, 0.8807, &
3015 0.8805, 0.8802, 0.8796, 0.8786, 0.8770, 0.8747, 0.8712, 0.8659, &
3016 0.8582, 0.8473, 0.8329, 0.8153, 0.7945, 0.7697, 0.7394, 0.7024, &
3017 0.6588, 0.6105, 0.8427, 0.8427, 0.8427, 0.8426, 0.8424, 0.8422, &
3018 0.8417, 0.8409, 0.8397, 0.8378, 0.8350, 0.8306, 0.8241, 0.8148, &
3019 0.8023, 0.7866, 0.7676, 0.7444, 0.7154, 0.6796, 0.6370, 0.5897, &
3020 0.7990, 0.7990, 0.7990, 0.7989, 0.7988, 0.7987, 0.7983, 0.7978, &
3021 0.7969, 0.7955, 0.7933, 0.7899, 0.7846, 0.7769, 0.7664, 0.7528, &
3022 0.7357, 0.7141, 0.6866, 0.6520, 0.6108, 0.5646/
3023
3024 data ((cah(i,j),i=1,22),j=11,15)/ &
3025 0.7515, 0.7515, 0.7515, 0.7515, 0.7514, 0.7513, 0.7511, 0.7507, &
3026 0.7501, 0.7491, 0.7476, 0.7450, 0.7409, 0.7347, 0.7261, 0.7144, &
3027 0.6992, 0.6793, 0.6533, 0.6203, 0.5805, 0.5357, 0.7020, 0.7020, &
3028 0.7020, 0.7019, 0.7019, 0.7018, 0.7017, 0.7015, 0.7011, 0.7005, &
3029 0.6993, 0.6974, 0.6943, 0.6894, 0.6823, 0.6723, 0.6588, 0.6406, &
3030 0.6161, 0.5847, 0.5466, 0.5034, 0.6518, 0.6518, 0.6518, 0.6518, &
3031 0.6518, 0.6517, 0.6517, 0.6515, 0.6513, 0.6508, 0.6500, 0.6485, &
3032 0.6459, 0.6419, 0.6359, 0.6273, 0.6151, 0.5983, 0.5755, 0.5458, &
3033 0.5095, 0.4681, 0.6017, 0.6017, 0.6017, 0.6017, 0.6016, 0.6016, &
3034 0.6016, 0.6015, 0.6013, 0.6009, 0.6002, 0.5989, 0.5967, 0.5932, &
3035 0.5879, 0.5801, 0.5691, 0.5535, 0.5322, 0.5043, 0.4700, 0.4308, &
3036 0.5518, 0.5518, 0.5518, 0.5518, 0.5518, 0.5518, 0.5517, 0.5516, &
3037 0.5514, 0.5511, 0.5505, 0.5493, 0.5473, 0.5441, 0.5393, 0.5322, &
3038 0.5220, 0.5076, 0.4878, 0.4617, 0.4297, 0.3929/
3039
3040 data ((cah(i,j),i=1,22),j=16,19)/ &
3041 0.5031, 0.5031, 0.5031, 0.5031, 0.5031, 0.5030, 0.5030, 0.5029, &
3042 0.5028, 0.5025, 0.5019, 0.5008, 0.4990, 0.4960, 0.4916, 0.4850, &
3043 0.4757, 0.4624, 0.4441, 0.4201, 0.3904, 0.3564, 0.4565, 0.4565, &
3044 0.4565, 0.4564, 0.4564, 0.4564, 0.4564, 0.4563, 0.4562, 0.4559, &
3045 0.4553, 0.4544, 0.4527, 0.4500, 0.4460, 0.4400, 0.4315, 0.4194, &
3046 0.4028, 0.3809, 0.3538, 0.3227, 0.4122, 0.4122, 0.4122, 0.4122, &
3047 0.4122, 0.4122, 0.4122, 0.4121, 0.4120, 0.4117, 0.4112, 0.4104, &
3048 0.4089, 0.4065, 0.4029, 0.3976, 0.3900, 0.3792, 0.3643, 0.3447, &
3049 0.3203, 0.2923, 0.3696, 0.3696, 0.3696, 0.3696, 0.3696, 0.3696, &
3050 0.3695, 0.3695, 0.3694, 0.3691, 0.3687, 0.3680, 0.3667, 0.3647, &
3051 0.3615, 0.3570, 0.3504, 0.3409, 0.3279, 0.3106, 0.2892, 0.2642/
3052
3053 !********************************************************************
3054 !-----table look-up for the reduction of clear-sky solar
3055 ! radiation due to co2. The fraction 0.0343 is the
3056 ! extraterrestrial solar flux in the co2 bands.
3057
3058 do k= 2, np+1
3059 do j= 1, n
3060 do i= 1, m
3061 xx=1./.3
3062 clog=log10(swc(i,j,k)*csm(i,j))
3063 wlog=log10(swh(i,j,k)*csm(i,j))
3064 ic=int( (clog+3.15)*xx+1.)
3065 iw=int( (wlog+4.15)*xx+1.)
3066 if(ic.lt.2)ic=2
3067 if(iw.lt.2)iw=2
3068 if(ic.gt.22)ic=22
3069 if(iw.gt.19)iw=19
3070 dc=clog-float(ic-2)*.3+3.
3071 dw=wlog-float(iw-2)*.3+4.
3072 x1=cah(1,iw-1)+(cah(1,iw)-cah(1,iw-1))*xx*dw
3073 x2=cah(ic-1,iw-1)+(cah(ic-1,iw)-cah(ic-1,iw-1))*xx*dw
3074 y2=x2+(cah(ic,iw-1)-cah(ic-1,iw-1))*xx*dc
3075 if (x1.lt.y2) x1=y2
3076 df(i,j,k)=df(i,j,k)+0.0343*(x1-y2)
3077 enddo
3078 enddo
3079 enddo
3080
3081 end subroutine flxco2
3082
3083 !*****************************************************************
3084
3085 subroutine o3prof (np, pres, ozone, its, ite, kts, kte, p, o3)
3086
3087 !*****************************************************************
3088 implicit none
3089 !*****************************************************************
3090 !
3091 integer iprof,m,np,its,ite,kts,kte
3092 integer i,k,ko,kk
3093 real pres(np),ozone(np)
3094 real p(its:ite,kts:kte),o3(its:ite,kts:kte)
3095
3096 ! Statement function
3097
3098 real Linear, x1, y1, x2, y2, x
3099 Linear(x1, y1, x2, y2, x) = &
3100 (y1 * (x2 - x) + y2 * (x - x1)) / (x2 - x1)
3101 !
3102 do k = 1,np
3103 pres(k) = alog(pres(k))
3104 enddo
3105 do k = kts,kte
3106 do i = its, ite
3107 p(i,k) = alog(p(i,k))
3108 end do
3109 end do
3110
3111 ! assume the pressure at model top is greater than pres(1)
3112 ! if it is not, this part needs to change
3113
3114 do i = its, ite
3115 ko = 1
3116 do k = kts+1, kte
3117 do while (ko .lt. np .and. p(i,k) .gt. pres(ko))
3118 ko = ko + 1
3119 end do
3120 o3(i,k) = Linear (pres(ko), ozone(ko), &
3121 pres(ko-1), ozone(ko-1), &
3122 p(i,k))
3123 ko = ko - 1
3124 end do
3125 end do
3126
3127 ! calculate top lay O3
3128
3129 do i = its, ite
3130 ko = 1
3131 k = kts
3132 do while (ko .le. np .and. p(i,k) .gt. pres(ko))
3133 ko = ko + 1
3134 end do
3135 IF (ko-1 .le. 1) then
3136 O3(i,k)=ozone(k)
3137 ELSE
3138 O3(i,k)=0.
3139 do kk=ko-2,1,-1
3140 O3(i,k)=O3(i,k)+ozone(kk)*(pres(kk+1)-pres(kk))
3141 enddo
3142 O3(i,k)=O3(i,k)/(pres(ko-1)-pres(1))
3143 ENDIF
3144 ! print*,'O3=',i,k,ko,O3(i,k),p(i,k),ko,pres(ko),pres(ko-1)
3145 end do
3146
3147 end subroutine o3prof
3148
3149 !-----------------------------------------
3150 SUBROUTINE gsfc_swinit(cen_lat, allowed_to_read)
3151
3152 REAL, INTENT(IN ) :: cen_lat
3153 LOGICAL, INTENT(IN ) :: allowed_to_read
3154
3155 center_lat=cen_lat
3156
3157 END SUBROUTINE gsfc_swinit
3158
3159
3160 END MODULE module_ra_gsfcsw