/[MITgcm]/MITgcm/pkg/fizhi/fizhi_swrad.F
ViewVC logotype

Contents of /MITgcm/pkg/fizhi/fizhi_swrad.F

Parent Directory Parent Directory | Revision Log Revision Log | View Revision Graph Revision Graph


Revision 1.25 - (show annotations) (download)
Fri Jun 17 16:51:24 2005 UTC (19 years ago) by ce107
Branch: MAIN
CVS Tags: checkpoint58l_post, checkpoint57t_post, checkpoint57o_post, checkpoint58e_post, checkpoint57v_post, checkpoint58u_post, checkpoint58w_post, checkpoint57m_post, checkpoint57s_post, checkpoint57k_post, checkpoint60, checkpoint61, checkpoint58r_post, checkpoint57i_post, checkpoint57y_post, checkpoint58n_post, checkpoint58x_post, checkpoint58t_post, checkpoint58h_post, checkpoint57y_pre, checkpoint58q_post, checkpoint59q, checkpoint59p, checkpoint59r, checkpoint58j_post, checkpoint59e, checkpoint59d, checkpoint59g, checkpoint59f, checkpoint59a, checkpoint59c, checkpoint59b, checkpoint59m, checkpoint59l, checkpoint59o, checkpoint59n, checkpoint59i, checkpoint59h, checkpoint59k, checkpoint59j, checkpoint57r_post, checkpoint59, checkpoint58, checkpoint58f_post, checkpoint57x_post, checkpoint57n_post, checkpoint58d_post, checkpoint58c_post, checkpoint57w_post, checkpoint57p_post, checkpint57u_post, checkpoint58a_post, checkpoint58i_post, checkpoint57q_post, checkpoint58g_post, checkpoint58o_post, checkpoint57z_post, checkpoint58y_post, checkpoint58k_post, checkpoint58v_post, checkpoint58s_post, checkpoint61f, checkpoint61g, checkpoint61d, checkpoint61e, checkpoint61b, checkpoint61c, checkpoint58p_post, checkpoint61a, checkpoint61j, checkpoint61k, checkpoint61h, checkpoint61i, checkpoint57j_post, checkpoint58b_post, checkpoint58m_post, checkpoint57l_post
Changes since 1.24: +12 -8 lines
Fixed length of lines with _d 0 additions to 72 chars

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/fizhi_swrad.F,v 1.24 2005/06/16 16:46:12 ce107 Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 subroutine swrio (nymd,nhms,bi,bj,ndswr,myid,istrip,npcs,
6 . low_level,mid_level,im,jm,lm,
7 . pz,plz,plze,dpres,pkht,pkz,tz,qz,oz,co2,
8 . albvisdr,albvisdf,albnirdr,albnirdf,
9 . dtradsw,dtswclr,radswg,swgclr,
10 . fdifpar,fdirpar,osr,osrclr,
11 . ptop,nswcld,cldsw,cswmo,nswlz,swlz,
12 . lpnt,imstturb,qliqave,fccave,landtype,xlats,xlons)
13
14 implicit none
15
16 c Input Variables
17 c ---------------
18 integer nymd,nhms,bi,bj,ndswr,myid,istrip,npcs
19 integer mid_level,low_level
20 integer im,jm,lm
21 _RL ptop
22 _RL pz(im,jm),plz(im,jm,lm),plze(im,jm,lm+1),dpres(im,jm,lm)
23 _RL pkht(im,jm,lm+1),pkz(im,jm,lm)
24 _RL tz(im,jm,lm),qz(im,jm,lm)
25 _RL oz(im,jm,lm)
26 _RL co2
27 _RL albvisdr(im,jm),albvisdf(im,jm),albnirdr(im,jm)
28 _RL albnirdf(im,jm)
29 _RL radswg(im,jm),swgclr(im,jm),fdifpar(im,jm),fdirpar(im,jm)
30 _RL osr(im,jm),osrclr(im,jm),dtradsw(im,jm,lm)
31 _RL dtswclr(im,jm,lm)
32 integer nswcld,nswlz
33 _RL cldsw(im,jm,lm),cswmo(im,jm,lm),swlz(im,jm,lm)
34 logical lpnt
35 integer imstturb
36 _RL qliqave(im,jm,lm),fccave(im,jm,lm)
37 integer landtype(im,jm)
38 _RL xlats(im,jm),xlons(im,jm)
39
40 c Local Variables
41 c ---------------
42 integer i,j,L,nn,nsecf
43 integer ntmstp,nymd2,nhms2
44 _RL getcon,grav,cp,undef
45 _RL ra,alf,reffw,reffi,tminv
46
47 parameter ( reffw = 10.0 )
48 parameter ( reffi = 65.0 )
49
50 _RL tdry(im,jm,lm)
51 _RL TEMP1(im,jm)
52 _RL TEMP2(im,jm)
53 _RL zenith (im,jm)
54 _RL cldtot (im,jm,lm)
55 _RL cldmxo (im,jm,lm)
56 _RL totcld (im,jm)
57 _RL cldlow (im,jm)
58 _RL cldmid (im,jm)
59 _RL cldhi (im,jm)
60 _RL taulow (im,jm)
61 _RL taumid (im,jm)
62 _RL tauhi (im,jm)
63 _RL tautype(im,jm,lm,3)
64 _RL tau(im,jm,lm)
65 _RL albedo(im,jm)
66
67 _RL PK(ISTRIP,lm)
68 _RL qzl(ISTRIP,lm),CLRO(ISTRIP,lm)
69 _RL TZL(ISTRIP,lm)
70 _RL OZL(ISTRIP,lm)
71 _RL PLE(ISTRIP,lm+1)
72 _RL COSZ(ISTRIP)
73 _RL dpstrip(ISTRIP,lm)
74
75 _RL albuvdr(ISTRIP),albuvdf(ISTRIP)
76 _RL albirdr(ISTRIP),albirdf(ISTRIP)
77 _RL difpar (ISTRIP),dirpar (ISTRIP)
78
79 _RL fdirir(istrip),fdifir(istrip)
80 _RL fdiruv(istrip),fdifuv(istrip)
81
82 _RL flux(istrip,lm+1)
83 _RL fluxclr(istrip,lm+1)
84 _RL dtsw(istrip,lm)
85 _RL dtswc(istrip,lm)
86
87 _RL taul (istrip,lm)
88 _RL reff (istrip,lm,2)
89 _RL tauc (istrip,lm,2)
90 _RL taua (istrip,lm)
91 _RL tstrip (istrip)
92 #ifdef ALLOW_DIAGNOSTICS
93 logical diagnostics_is_on
94 external diagnostics_is_on
95 _RL tmpdiag(im,jm),tmpdiag2(im,jm)
96 #endif
97
98 logical first
99 data first /.true./
100
101 C **********************************************************************
102 C **** INITIALIZATION ****
103 C **********************************************************************
104
105 grav = getcon('GRAVITY')
106 cp = getcon('CP')
107 undef = getcon('UNDEF')
108
109 NTMSTP = nsecf(NDSWR)
110 TMINV = 1./float(ntmstp)
111
112 C Compute Temperature from Theta
113 C ------------------------------
114 do L=1,lm
115 do j=1,jm
116 do i=1,im
117 tdry(i,j,L) = tz(i,j,L)*pkz(i,j,L)
118 enddo
119 enddo
120 enddo
121
122 if (first .and. myid.eq.1 ) then
123 print *
124 print *,'Low-Level Clouds are Grouped between levels: ',
125 . lm,' and ',low_level
126 print *,'Mid-Level Clouds are Grouped between levels: ',
127 . low_level-1,' and ',mid_level
128 print *
129 first = .false.
130 endif
131
132 C **********************************************************************
133 C **** CALCULATE COSINE OF THE ZENITH ANGLE ****
134 C **********************************************************************
135
136 CALL ASTRO ( NYMD, NHMS, XLATS,XLONS, im*jm, TEMP1,RA )
137 NYMD2 = NYMD
138 NHMS2 = NHMS
139 CALL TICK ( NYMD2, NHMS2, NTMSTP )
140 CALL ASTRO ( NYMD2, NHMS2, XLATS,XLONS, im*jm, TEMP2,RA )
141
142 do j = 1,jm
143 do i = 1,im
144 zenith(I,j) = TEMP1(I,j) + TEMP2(I,j)
145 IF( zenith(I,j) .GT. 0.0 ) THEN
146 zenith(I,j) = (2./3.)*( zenith(I,j)-TEMP2(I,j)*TEMP1(I,j)
147 . / zenith(I,j) )
148 ENDIF
149 ENDDO
150 ENDDO
151
152 C **********************************************************************
153 c **** Compute Two-Dimension Total Cloud Fraction (0-1) ****
154 C **********************************************************************
155
156 c Initialize Clear Sky Probability for Low, Mid, and High Clouds
157 c --------------------------------------------------------------
158 do j =1,jm
159 do i =1,im
160 cldlow(i,j) = 0.0
161 cldmid(i,j) = 0.0
162 cldhi (i,j) = 0.0
163 enddo
164 enddo
165
166 c Adjust cloud fractions and cloud liquid water due to moist turbulence
167 c ---------------------------------------------------------------------
168 if(imstturb.ne.0) then
169 do L =1,lm
170 do j =1,jm
171 do i =1,im
172 cldtot(i,j,L)=min(1.0 _d 0,max(cldsw(i,j,L),fccave(i,j,L)/
173 $ imstturb))
174 cldmxo(i,j,L)=min(1.0 _d 0,cswmo(i,j,L))
175 swlz(i,j,L)=swlz(i,j,L)+qliqave(i,j,L)/imstturb
176 enddo
177 enddo
178 enddo
179 else
180 do L =1,lm
181 do j =1,jm
182 do i =1,im
183 cldtot(i,j,L) = min( 1.0 _d 0,cldsw(i,j,L) )
184 cldmxo(i,j,L) = min( 1.0 _d 0,cswmo(i,j,L) )
185 enddo
186 enddo
187 enddo
188 endif
189
190 c Compute 3-D Cloud Fractions
191 c ---------------------------
192 if( nswcld.ne.0 ) then
193 do L = 1,lm
194 do j = 1,jm
195 do i = 1,im
196 c Compute Low-Mid-High Maximum Overlap Cloud Fractions
197 c ----------------------------------------------------
198 if( L.lt.mid_level ) then
199 cldhi (i,j) = max( cldtot(i,j,L),cldhi (i,j) )
200 elseif( L.ge.mid_level .and.
201 . L.lt.low_level ) then
202 cldmid(i,j) = max( cldtot(i,j,L),cldmid(i,j) )
203 elseif( L.ge.low_level ) then
204 cldlow(i,j) = max( cldtot(i,j,L),cldlow(i,j) )
205 endif
206
207 enddo
208 enddo
209 enddo
210 endif
211
212 c Totcld => Product of Clear Sky Probabilities for Low, Mid, and High Clouds
213 c --------------------------------------------------------------------------
214 do j = 1,jm
215 do i = 1,im
216 totcld(i,j) = 1.0 - (1.-cldhi (i,j))
217 . * (1.-cldmid(i,j))
218 . * (1.-cldlow(i,j))
219 enddo
220 enddo
221
222 #ifdef ALLOW_DIAGNOSTICS
223
224 c Compute Cloud Diagnostics
225 c -------------------------
226 if(diagnostics_is_on('CLDFRC ',myid) ) then
227 call diagnostics_fill(totcld,'CLDFRC ',0,1,3,bi,bj,myid)
228 endif
229
230 if(diagnostics_is_on('CLDRAS ',myid) ) then
231 do L=1,lm
232 do j=1,jm
233 do i=1,im
234 tmpdiag(i,j) = cswmo(i,j,L)
235 enddo
236 enddo
237 call diagnostics_fill(tmpdiag,'CLDRAS ',L,1,3,bi,bj,myid)
238 enddo
239 endif
240
241 if(diagnostics_is_on('CLDTOT ',myid) ) then
242 do L=1,lm
243 do j=1,jm
244 do i=1,im
245 tmpdiag(i,j) = cldtot(i,j,L)
246 enddo
247 enddo
248 call diagnostics_fill(tmpdiag,'CLDTOT ',L,1,3,bi,bj,myid)
249 enddo
250 endif
251
252 if(diagnostics_is_on('CLDLOW ',myid) ) then
253 call diagnostics_fill(cldlow,'CLDLOW ',0,1,3,bi,bj,myid)
254 endif
255
256 if(diagnostics_is_on('CLDMID ',myid) ) then
257 call diagnostics_fill(cldmid,'CLDMID ',0,1,3,bi,bj,myid)
258 endif
259
260 if(diagnostics_is_on('CLDHI ',myid) ) then
261 call diagnostics_fill(cldhi,'CLDHI ',0,1,3,bi,bj,myid)
262 endif
263
264 if(diagnostics_is_on('LZRAD ',myid) ) then
265 do L=1,lm
266 do j=1,jm
267 do i=1,im
268 tmpdiag(i,j) = swlz(i,j,L) * 1.0e6
269 enddo
270 enddo
271 call diagnostics_fill(tmpdiag,'LZRAD ',L,1,3,bi,bj,myid)
272 enddo
273 endif
274
275 c Albedo Diagnostics
276 c ------------------
277 if(diagnostics_is_on('ALBVISDR',myid) ) then
278 call diagnostics_fill(albvisdr,'ALBVISDR',0,1,3,bi,bj,myid)
279 endif
280
281 if(diagnostics_is_on('ALBVISDF',myid) ) then
282 call diagnostics_fill(albvisdf,'ALBVISDF',0,1,3,bi,bj,myid)
283 endif
284
285 if(diagnostics_is_on('ALBNIRDR',myid) ) then
286 call diagnostics_fill(albnirdr,'ALBNIRDR',0,1,3,bi,bj,myid)
287 endif
288
289 if(diagnostics_is_on('ALBNIRDF',myid) ) then
290 call diagnostics_fill(albnirdf,'ALBNIRDF',0,1,3,bi,bj,myid)
291 endif
292
293 #endif
294
295 C Compute Optical Thicknesses and Diagnostics
296 C -------------------------------------------
297 call opthk(tdry,plz,plze,swlz,cldtot,cldmxo,landtype,im,jm,lm,
298 . tautype)
299
300 do L = 1,lm
301 do j = 1,jm
302 do i = 1,im
303 tau(i,j,L) = tautype(i,j,L,1)+tautype(i,j,L,2)+tautype(i,j,L,3)
304 enddo
305 enddo
306 enddo
307
308 #ifdef ALLOW_DIAGNOSTICS
309 if(diagnostics_is_on('TAUAVE ',myid) ) then
310 do L=1,lm
311 do j=1,jm
312 do i=1,im
313 tmpdiag(i,j) = tau(i,j,L)*100/(plze(i,j,L+1)-plze(i,j,L))
314 enddo
315 enddo
316 call diagnostics_fill(tmpdiag,'TAUAVE ',L,1,3,bi,bj,myid)
317 enddo
318 endif
319
320 if(diagnostics_is_on('TAUCLD ',myid) .and.
321 . diagnostics_is_on('TAUCLDC ',myid) ) then
322 do L=1,lm
323 do j=1,jm
324 do i=1,im
325 if( cldtot(i,j,L).ne.0.0 ) then
326 tmpdiag(i,j) = tau(i,j,L)*100/(plze(i,j,L+1)-plze(i,j,L))
327 tmpdiag2(i,j) = 1.
328 else
329 tmpdiag(i,j) = 0.
330 tmpdiag2(i,j) = 0.
331 endif
332 enddo
333 enddo
334 call diagnostics_fill(tmpdiag,'TAUCLD ',L,1,3,bi,bj,myid)
335 call diagnostics_fill(tmpdiag,'TAUCLDC ',L,1,3,bi,bj,myid)
336 enddo
337 endif
338
339 c Compute Low, Mid, and High Cloud Optical Depth Diagnostics
340 c ----------------------------------------------------------
341 if(diagnostics_is_on('TAULOW ',myid) .and.
342 . diagnostics_is_on('TAULOWC ',myid) ) then
343 do j = 1,jm
344 do i = 1,im
345 taulow(i,j) = 0.0
346 if( cldlow(i,j).ne.0.0 ) then
347 do L = low_level,lm
348 taulow(i,j) = taulow(i,j) + tau(i,j,L)
349 enddo
350 tmpdiag2(i,j) = 1.
351 else
352 tmpdiag(i,j) = 0.
353 endif
354 enddo
355 enddo
356 call diagnostics_fill(taulow,'TAULOW ',0,1,3,bi,bj,myid)
357 call diagnostics_fill(tmpdiag2,'TAULOWC ',0,1,3,bi,bj,myid)
358 endif
359
360 if(diagnostics_is_on('TAUMID ',myid) .and.
361 . diagnostics_is_on('TAUMIDC ',myid) ) then
362 do j = 1,jm
363 do i = 1,im
364 taumid(i,j) = 0.0
365 if( cldmid(i,j).ne.0.0 ) then
366 do L = mid_level,low_level+1
367 taumid(i,j) = taumid(i,j) + tau(i,j,L)
368 enddo
369 tmpdiag2(i,j) = 1.
370 else
371 tmpdiag(i,j) = 0.
372 endif
373 enddo
374 enddo
375 call diagnostics_fill(taumid,'TAUMID ',0,1,3,bi,bj,myid)
376 call diagnostics_fill(tmpdiag2,'TAUMIDC ',0,1,3,bi,bj,myid)
377 endif
378
379 if(diagnostics_is_on('TAUHI ',myid) .and.
380 . diagnostics_is_on('TAUHIC ',myid) ) then
381 do j = 1,jm
382 do i = 1,im
383 tauhi(i,j) = 0.0
384 if( cldhi(i,j).ne.0.0 ) then
385 do L = 1,mid_level+1
386 tauhi(i,j) = tauhi(i,j) + tau(i,j,L)
387 enddo
388 tmpdiag2(i,j) = 1.
389 else
390 tmpdiag(i,j) = 0.
391 endif
392 enddo
393 enddo
394 call diagnostics_fill(tauhi,'TAUHI ',0,1,3,bi,bj,myid)
395 call diagnostics_fill(tmpdiag2,'TAUHIC ',0,1,3,bi,bj,myid)
396 endif
397 #endif
398
399 C***********************************************************************
400 C **** LOOP OVER GLOBAL REGIONS ****
401 C **********************************************************************
402
403 do 1000 nn = 1,npcs
404
405 C **********************************************************************
406 C **** VARIABLE INITIALIZATION ****
407 C **********************************************************************
408
409 CALL STRIP ( zenith,COSZ,im*jm,ISTRIP,1,NN )
410
411 CALL STRIP ( plze, ple ,im*jm,ISTRIP,lm+1,NN)
412 CALL STRIP ( pkz , pk ,im*jm,ISTRIP,lm ,NN)
413 CALL STRIP ( dpres,dpstrip,im*jm,ISTRIP,lm ,NN)
414 CALL STRIP ( tdry, tzl ,im*jm,ISTRIP,lm ,NN)
415 CALL STRIP ( qz , qzl ,im*jm,ISTRIP,lm ,NN)
416 CALL STRIP ( oz , ozl ,im*jm,ISTRIP,lm ,NN)
417 CALL STRIP ( tau , taul ,im*jm,ISTRIP,lm ,NN)
418
419 CALL STRIP ( albvisdr,albuvdr,im*jm,ISTRIP,1,NN )
420 CALL STRIP ( albvisdf,albuvdf,im*jm,ISTRIP,1,NN )
421 CALL STRIP ( albnirdr,albirdr,im*jm,ISTRIP,1,NN )
422 CALL STRIP ( albnirdf,albirdf,im*jm,ISTRIP,1,NN )
423
424 call strip ( cldtot,clro,im*jm,istrip,lm,nn )
425
426 c Partition Tau between Water and Ice Particles
427 c ---------------------------------------------
428 do L= 1,lm
429 do i= 1,istrip
430 alf = min( max((tzl(i,l)-253.15)/20.,0. _d 0) ,1. _d 0)
431 taua(i,L) = 0.
432
433 if( alf.ne.0.0 .and. alf.ne.1.0 ) then
434 tauc(i,L,1) = taul(i,L)/(1.+alf/(1-alf) * (reffi/reffw*0.80) )
435 tauc(i,L,2) = taul(i,L)/(1.+(1-alf)/alf * (reffw/reffi*1.25) )
436 endif
437
438 if( alf.eq.0.0 ) then
439 tauc(i,L,1) = taul(i,L)
440 tauc(i,L,2) = 0.0
441 endif
442
443 if( alf.eq.1.0 ) then
444 tauc(i,L,1) = 0.0
445 tauc(i,L,2) = taul(i,L)
446 endif
447
448 reff(i,L,1) = reffi
449 reff(i,L,2) = reffw
450 enddo
451 enddo
452
453 call sorad ( istrip,1,1,lm,ple,tzl,qzl,ozl,co2,
454 . tauc,reff,clro,mid_level,low_level,taua,
455 . albirdr,albirdf,albuvdr,albuvdf,cosz,
456 . flux,fluxclr,fdirir,fdifir,dirpar,difpar,
457 . fdiruv,fdifuv )
458
459 C **********************************************************************
460 C **** Compute Mass-Weighted Theta Tendencies from Fluxes ****
461 C **********************************************************************
462
463 do l=1,lm
464 do i=1,istrip
465 alf = grav*(ple(i,Lm+1)-ptop)/(cp*dpstrip(i,L)*100)
466 dtsw (i,L) = alf*( flux (i,L)-flux (i,L+1) )/pk(i,L)
467 dtswc(i,L) = alf*( fluxclr(i,L)-fluxclr(i,L+1) )/pk(i,L)
468 enddo
469 enddo
470
471 call paste ( dtsw , dtradsw ,istrip,im*jm,lm,nn )
472 call paste ( dtswc, dtswclr ,istrip,im*jm,lm,nn )
473
474 call paste ( flux (1,1),osr ,istrip,im*jm,1,nn )
475 call paste ( fluxclr(1,1),osrclr,istrip,im*jm,1,nn )
476
477 call paste ( flux (1,lm+1),radswg,istrip,im*jm,1,nn )
478 call paste ( fluxclr(1,lm+1),swgclr,istrip,im*jm,1,nn )
479
480 call paste ( dirpar,fdirpar,istrip,im*jm,1,nn )
481 call paste ( difpar,fdifpar,istrip,im*jm,1,nn )
482
483 c Calculate Mean Albedo
484 c ---------------------
485 do i=1,istrip
486 if( cosz(i).gt.0.0 ) then
487 tstrip(i) = 1.0 - flux(i,lm+1)/
488 . ( fdirir(i)+fdifir(i)+dirpar(i)+difpar(i) + fdiruv(i)+fdifuv(i) )
489 if( tstrip(i).lt.0.0 ) tstrip(i) = undef
490 else
491 tstrip(i) = undef
492 endif
493 enddo
494 call paste ( tstrip,albedo,istrip,im*jm,1,nn )
495
496 1000 CONTINUE
497
498 #ifdef ALLOW_DIAGNOSTICS
499
500 c Mean Albedo Diagnostic
501 c ----------------------
502 if(diagnostics_is_on('ALBEDO ',myid) .and.
503 . diagnostics_is_on('ALBEDOC ',myid) ) then
504 do j=1,jm
505 do i=1,im
506 if( albedo(i,j).ne.undef ) then
507 tmpdiag(i,j) = albedo(i,j)
508 tmpdiag2(i,j) = 1.
509 else
510 tmpdiag(i,j) = 0.
511 tmpdiag2(i,j) = 0.
512 endif
513 enddo
514 enddo
515 call diagnostics_fill(tmpdiag,'ALBEDO ',0,1,3,bi,bj,myid)
516 call diagnostics_fill(tmpdiag2,'ALBEDOC ',0,1,3,bi,bj,myid)
517 endif
518 #endif
519
520 C **********************************************************************
521 C **** ZERO-OUT OR BUMP COUNTERS ****
522 C **********************************************************************
523
524 nswlz = 0
525 nswcld = 0
526 imstturb = 0
527
528 do L = 1,lm
529 do j = 1,jm
530 do i = 1,im
531 fccave(i,j,L) = 0.
532 qliqave(i,j,L) = 0.
533 enddo
534 enddo
535 enddo
536
537 return
538 end
539 subroutine opthk ( tl,pl,ple,lz,cf,cfm,lwi,im,jm,lm,tau )
540 C***********************************************************************
541 C
542 C PURPOSE:
543 C ========
544 C Compute cloud optical thickness using temperature and
545 C cloud fractions.
546 C
547 C INPUT:
548 C ======
549 C tl ......... Temperature at Model Layers (K)
550 C pl ......... Model Layer Pressures (mb)
551 C ple ........ Model Edge Pressures (mb)
552 C lz ......... Diagnosed Convective and Large-Scale Cloud Water (g/g)
553 C cf ......... Total Cloud Fraction (Random + Maximum Overlap)
554 C cfm ........ Maximum Overlap Cloud Fraction
555 C lwi ........ Surface Land Type
556 C im ......... Longitudinal dimension
557 C jm ......... Latitudinal dimension
558 C lm ......... Vertical dimension
559 C
560 C OUTPUT:
561 C =======
562 C tau ........ Cloud Optical Thickness (non-dimensional)
563 C tau(im,jm,lm,1): Suspended Ice
564 C tau(im,jm,lm,2): Suspended Water
565 C tau(im,jm,lm,3): Raindrops
566 C
567 C***********************************************************************
568
569 implicit none
570
571 integer im,jm,lm,i,j,L
572
573 _RL tl(im,jm,lm)
574 _RL pl(im,jm,lm)
575 _RL ple(im,jm,lm+1)
576 _RL lz(im,jm,lm)
577 _RL cf(im,jm,lm)
578 _RL cfm(im,jm,lm)
579 _RL tau(im,jm,lm,3)
580 integer lwi(im,jm)
581
582 _RL dp, alf, fracls, fraccu
583 _RL tauice, tauh2o, tauras
584
585 c Compute Cloud Optical Depths
586 c ----------------------------
587 do L=1,lm
588 do j=1,jm
589 do i=1,im
590 alf = min( max(( tl(i,j,L)-233.15)/20.,0. _d 0), 1. _d 0)
591 dp = ple(i,j,L+1)-ple(i,j,L)
592 tau(i,j,L,1) = 0.0
593 tau(i,j,L,2) = 0.0
594 tau(i,j,L,3) = 0.0
595 if( cf(i,j,L).ne.0.0 ) then
596
597 c Determine fraction of large-scale and cumulus clouds
598 c ----------------------------------------------------
599 fracls = ( cf(i,j,L)-cfm(i,j,L) )/cf(i,j,L)
600 fraccu = 1.0-fracls
601
602 c Define tau for large-scale ice and water clouds
603 c Note: tauice is scaled between (0.02 & .2) for: 0 < lz < 2 mg/kg over Land
604 c Note: tauh2o is scaled between (0.20 & 20) for: 0 < lz < 5 mg/kg over Land
605 c Note: tauh2o is scaled between (0.20 & 12) for: 0 < lz < 50 mg/kg over Ocean
606 c -------------------------------------------------------------------------------
607
608 c Large-Scale Ice
609 c ---------------
610 tauice = max( 0.0002 _d 0, 0.002*min(500*lz(i,j,L)*1000,
611 $ 1.0 _d 0) )
612 tau(i,j,L,1) = fracls*(1-alf)*tauice*dp
613
614 c Large-Scale Water
615 c -----------------
616 C Over Land
617 if( lwi(i,j).le.10 ) then
618 tauh2o = max( 0.0020 _d 0, 0.200*min(200*lz(i,j,L)*1000,
619 $ 1.0 _d 0) )
620 tau(i,j,L,3) = fracls*alf*tauh2o*dp
621 else
622 C Non-Precipitation Clouds Over Ocean
623 if( lz(i,j,L).eq.0.0 ) then
624 tauh2o = .12
625 tau(i,j,L,2) = fracls*alf*tauh2o*dp
626 else
627 C Over Ocean
628 tauh2o = max( 0.0020 _d 0, 0.120*min( 20*lz(i,j,L)*1000,
629 $ 1.0 _d 0) )
630 tau(i,j,L,3) = fracls*alf*tauh2o*dp
631 endif
632 endif
633
634 c Sub-Grid Convective
635 c -------------------
636 if( tl(i,j,L).gt.210.0 ) then
637 tauras = .16
638 tau(i,j,L,3) = tau(i,j,L,3) + fraccu*tauras*dp
639 else
640 tauras = tauice
641 tau(i,j,L,1) = tau(i,j,L,1) + fraccu*tauras*dp
642 endif
643
644 c Define tau for large-scale and cumulus clouds
645 c ---------------------------------------------
646 ccc tau(i,j,L) = ( fracls*( alf*tauh2o + (1-alf)*tauice )
647 ccc . + fraccu*tauras )*dp
648
649 endif
650
651 enddo
652 enddo
653 enddo
654
655 return
656 end
657 subroutine sorad(m,n,ndim,np,pl,ta,wa,oa,co2,
658 * taucld,reff,fcld,ict,icb,
659 * taual,rsirbm,rsirdf,rsuvbm,rsuvdf,cosz,
660 * flx,flc,fdirir,fdifir,fdirpar,fdifpar,
661 * fdiruv,fdifuv)
662 c********************************************************************
663 c
664 c This routine computes solar fluxes due to the absoption by water
665 c vapor, ozone, co2, o2, clouds, and aerosols and due to the
666 c scattering by clouds, aerosols, and gases.
667 c
668 c This is a vectorized code. It computes the fluxes simultaneous for
669 c (m x n) soundings, which is a subset of the (m x ndim) soundings.
670 c In a global climate model, m and ndim correspond to the numbers of
671 c grid boxes in the zonal and meridional directions, respectively.
672 c
673 c Ice and liquid cloud particles are allowed to co-exist in any of the
674 c np layers. Two sets of cloud parameters are required as inputs, one
675 c for ice paticles and the other for liquid particles. These parameters
676 c are optical thickness (taucld) and effective particle size (reff).
677 c
678 c If no information is available for reff, a default value of
679 c 10 micron for liquid water and 75 micron for ice can be used.
680 c
681 c Clouds are grouped into high, middle, and low clouds separated by the
682 c level indices ict and icb. For detail, see the subroutine cldscale.
683 c
684 c----- Input parameters:
685 c units size
686 c number of soundings in zonal direction (m) n/d 1
687 c number of soundings in meridional direction (n) n/d 1
688 c maximum number of soundings in n/d 1
689 c meridional direction (ndim)
690 c number of atmospheric layers (np) n/d 1
691 c level pressure (pl) mb m*ndim*(np+1)
692 c layer temperature (ta) k m*ndim*np
693 c layer specific humidity (wa) gm/gm m*ndim*np
694 c layer ozone concentration (oa) gm/gm m*ndim*np
695 c co2 mixing ratio by volumn (co2) parts/part 1
696 c cloud optical thickness (taucld) n/d m*ndim*np*2
697 c index 1 for ice particles
698 c index 2 for liquid drops
699 c effective cloud-particle size (reff) micrometer m*ndim*np*2
700 c index 1 for ice particles
701 c index 2 for liquid drops
702 c cloud amount (fcld) fraction m*ndim*np
703 c level index separating high and middle n/d 1
704 c clouds (ict)
705 c level index separating middle and low clouds n/d 1
706 c clouds (icb)
707 c aerosol optical thickness (taual) n/d m*ndim*np
708 c solar ir surface albedo for beam fraction m*ndim
709 c radiation (rsirbm)
710 c solar ir surface albedo for diffuse fraction m*ndim
711 c radiation (rsirdf)
712 c uv + par surface albedo for beam fraction m*ndim
713 c radiation (rsuvbm)
714 c uv + par surface albedo for diffuse fraction m*ndim
715 c radiation (rsuvdf)
716 c cosine of solar zenith angle (cosz) n/d m*ndim
717 c
718 c----- Output parameters
719 c
720 c all-sky flux (downward minus upward) (flx) fraction m*ndim*(np+1)
721 c clear-sky flux (downward minus upward) (flc) fraction m*ndim*(np+1)
722 c all-sky direct downward ir (0.7-10 micron)
723 c flux at the surface (fdirir) fraction m*ndim
724 c all-sky diffuse downward ir flux at
725 c the surface (fdifir) fraction m*ndim
726 c all-sky direct downward par (0.4-0.7 micron)
727 c flux at the surface (fdirpar) fraction m*ndim
728 c all-sky diffuse downward par flux at
729 c the surface (fdifpar) fraction m*ndim
730 *
731 c----- Notes:
732 c
733 c (1) The unit of flux is fraction of the incoming solar radiation
734 c at the top of the atmosphere. Therefore, fluxes should
735 c be equal to flux multiplied by the extra-terrestrial solar
736 c flux and the cosine of solar zenith angle.
737 c (2) Clouds and aerosols can be included in any layers by specifying
738 c fcld(i,j,k), taucld(i,j,k,*) and taual(i,j,k), k=1,np.
739 c For an atmosphere without clouds and aerosols,
740 c set fcld(i,j,k)=taucld(i,j,k,*)=taual(i,j,k)=0.0.
741 c (3) Aerosol single scattering albedos and asymmetry
742 c factors are specified in the subroutines solir and soluv.
743 c (4) pl(i,j,1) is the pressure at the top of the model, and
744 c pl(i,j,np+1) is the surface pressure.
745 c (5) the pressure levels ict and icb correspond approximately
746 c to 400 and 700 mb.
747 c
748 c**************************************************************************
749
750 implicit none
751
752 c-----Explicit Inline Directives
753
754 #ifdef CRAY
755 #ifdef f77
756 cfpp$ expand (expmn)
757 #endif
758 #endif
759 _RL expmn
760
761 c-----input parameters
762
763 integer m,n,ndim,np,ict,icb
764 _RL pl(m,ndim,np+1),ta(m,ndim,np),wa(m,ndim,np),oa(m,ndim,np)
765 _RL taucld(m,ndim,np,2),reff(m,ndim,np,2)
766 _RL fcld(m,ndim,np),taual(m,ndim,np)
767 _RL rsirbm(m,ndim),rsirdf(m,ndim),
768 * rsuvbm(m,ndim),rsuvdf(m,ndim),cosz(m,ndim),co2
769
770 c-----output parameters
771
772 _RL flx(m,ndim,np+1),flc(m,ndim,np+1)
773 _RL fdirir(m,ndim),fdifir(m,ndim)
774 _RL fdirpar(m,ndim),fdifpar(m,ndim)
775 _RL fdiruv(m,ndim),fdifuv(m,ndim)
776
777 c-----temporary array
778
779 integer i,j,k
780 _RL cc(m,n,3),tauclb(m,n,np),tauclf(m,n,np)
781 _RL dp(m,n,np),wh(m,n,np),oh(m,n,np),scal(m,n,np)
782 _RL swh(m,n,np+1),so2(m,n,np+1),df(m,n,np+1)
783 _RL sdf(m,n),sclr(m,n),csm(m,n),x
784
785 c-----------------------------------------------------------------
786
787 do j= 1, n
788 do i= 1, m
789
790 swh(i,j,1)=0.
791 so2(i,j,1)=0.
792
793 c-----csm is the effective secant of the solar zenith angle
794 c see equation (12) of Lacis and Hansen (1974, JAS)
795
796 csm(i,j)=35./sqrt(1224.*cosz(i,j)*cosz(i,j)+1.)
797
798 enddo
799 enddo
800
801
802 do k= 1, np
803 do j= 1, n
804 do i= 1, m
805
806 c-----compute layer thickness and pressure-scaling function.
807 c indices for the surface level and surface layer
808 c are np+1 and np, respectively.
809
810 dp(i,j,k)=pl(i,j,k+1)-pl(i,j,k)
811 scal(i,j,k)=dp(i,j,k)*(.5*(pl(i,j,k)+pl(i,j,k+1))/300.)**.8
812
813 c-----compute scaled water vapor amount, unit is g/cm**2
814
815 wh(i,j,k)=1.02*wa(i,j,k)*scal(i,j,k)*
816 * (1.+0.00135*(ta(i,j,k)-240.))
817 swh(i,j,k+1)=swh(i,j,k)+wh(i,j,k)
818
819 c-----compute ozone amount, unit is (cm-atm)stp.
820
821 oh(i,j,k)=1.02*oa(i,j,k)*dp(i,j,k)*466.7
822
823 enddo
824 enddo
825 enddo
826
827
828 c-----scale cloud optical thickness in each layer from taucld (with
829 c cloud amount fcld) to tauclb and tauclf (with cloud amount cc).
830 c tauclb is the scaled optical thickness for beam radiation and
831 c tauclf is for diffuse radiation.
832
833 call cldscale(m,n,ndim,np,cosz,fcld,taucld,ict,icb,
834 * cc,tauclb,tauclf)
835
836 c-----initialize fluxes for all-sky (flx), clear-sky (flc), and
837 c flux reduction (df)
838
839 do k=1, np+1
840 do j=1, n
841 do i=1, m
842 flx(i,j,k)=0.
843 flc(i,j,k)=0.
844 df(i,j,k)=0.
845 enddo
846 enddo
847 enddo
848
849 c-----compute solar ir fluxes
850
851 call solir (m,n,ndim,np,wh,taucld,tauclb,tauclf,reff,ict,icb
852 * ,fcld,cc,taual,csm,rsirbm,rsirdf,flx,flc,fdirir,fdifir)
853
854 c-----compute and update uv and par fluxes
855
856 call soluv (m,n,ndim,np,oh,dp,taucld,tauclb,tauclf,reff,ict,icb
857 * ,fcld,cc,taual,csm,rsuvbm,rsuvdf,flx,flc
858 * ,fdirpar,fdifpar,fdiruv,fdifuv)
859
860 c-----compute scaled amount of o2 (so2), unit is (cm-atm)stp.
861
862 do k= 1, np
863 do j= 1, n
864 do i= 1, m
865 so2(i,j,k+1)=so2(i,j,k)+165.22*scal(i,j,k)
866 enddo
867 enddo
868 enddo
869
870 c-----compute flux reduction due to oxygen following
871 c chou (J. climate, 1990). The fraction 0.0287 is the
872 c extraterrestrial solar flux in the o2 bands.
873
874 do k= 2, np+1
875 do j= 1, n
876 do i= 1, m
877 x=so2(i,j,k)*csm(i,j)
878 df(i,j,k)=df(i,j,k)+0.0287*(1.-expmn(-0.00027*sqrt(x)))
879 enddo
880 enddo
881 enddo
882
883 c-----compute scaled amounts for co2 (so2). unit is (cm-atm)stp.
884
885 do k= 1, np
886 do j= 1, n
887 do i= 1, m
888 so2(i,j,k+1)=so2(i,j,k)+co2*789.*scal(i,j,k)
889 enddo
890 enddo
891 enddo
892
893 c-----compute and update flux reduction due to co2 following
894 c chou (J. Climate, 1990)
895
896 call flxco2(m,n,np,so2,swh,csm,df)
897
898 c-----adjust for the effect of o2 cnd co2 on clear-sky fluxes.
899
900 do k= 2, np+1
901 do j= 1, n
902 do i= 1, m
903 flc(i,j,k)=flc(i,j,k)-df(i,j,k)
904 enddo
905 enddo
906 enddo
907
908 c-----adjust for the all-sky fluxes due to o2 and co2. It is
909 c assumed that o2 and co2 have no effects on solar radiation
910 c below clouds. clouds are assumed randomly overlapped.
911
912 do j=1,n
913 do i=1,m
914 sdf(i,j)=0.0
915 sclr(i,j)=1.0
916 enddo
917 enddo
918
919 do k=1,np
920 do j=1,n
921 do i=1,m
922
923 c-----sclr is the fraction of clear sky.
924 c sdf is the flux reduction below clouds.
925
926 if(fcld(i,j,k).gt.0.01) then
927 sdf(i,j)=sdf(i,j)+df(i,j,k)*sclr(i,j)*fcld(i,j,k)
928 sclr(i,j)=sclr(i,j)*(1.-fcld(i,j,k))
929 endif
930 flx(i,j,k+1)=flx(i,j,k+1)-sdf(i,j)-df(i,j,k+1)*sclr(i,j)
931
932 enddo
933 enddo
934 enddo
935
936 c-----adjust for the direct downward ir flux.
937 do j= 1, n
938 do i= 1, m
939 x = (1.-rsirdf(i,j))*fdifir(i,j) + (1.-rsirbm(i,j))*fdirir(i,j)
940 x = (-sdf(i,j)-df(i,j,np+1)*sclr(i,j))/x
941 fdirir(i,j)=fdirir(i,j)*(1.+x)
942 fdifir(i,j)=fdifir(i,j)*(1.+x)
943 enddo
944 enddo
945
946 return
947 end
948
949 c********************************************************************
950
951 subroutine cldscale (m,n,ndim,np,cosz,fcld,taucld,ict,icb,
952 * cc,tauclb,tauclf)
953
954 c********************************************************************
955 c
956 c This subroutine computes the covers of high, middle, and
957 c low clouds and scales the cloud optical thickness.
958 c
959 c To simplify calculations in a cloudy atmosphere, clouds are
960 c grouped into high, middle and low clouds separated by the levels
961 c ict and icb (level 1 is the top of the atmosphere).
962 c
963 c Within each of the three groups, clouds are assumed maximally
964 c overlapped, and the cloud cover (cc) of a group is the maximum
965 c cloud cover of all the layers in the group. The optical thickness
966 c (taucld) of a given layer is then scaled to new values (tauclb and
967 c tauclf) so that the layer reflectance corresponding to the cloud
968 c cover cc is the same as the original reflectance with optical
969 c thickness taucld and cloud cover fcld.
970 c
971 c---input parameters
972 c
973 c number of grid intervals in zonal direction (m)
974 c number of grid intervals in meridional direction (n)
975 c maximum number of grid intervals in meridional direction (ndim)
976 c number of atmospheric layers (np)
977 c cosine of the solar zenith angle (cosz)
978 c fractional cloud cover (fcld)
979 c cloud optical thickness (taucld)
980 c index separating high and middle clouds (ict)
981 c index separating middle and low clouds (icb)
982 c
983 c---output parameters
984 c
985 c fractional cover of high, middle, and low clouds (cc)
986 c scaled cloud optical thickness for beam radiation (tauclb)
987 c scaled cloud optical thickness for diffuse radiation (tauclf)
988 c
989 c********************************************************************
990
991 implicit none
992
993 c-----input parameters
994
995 integer m,n,ndim,np,ict,icb
996 _RL cosz(m,ndim),fcld(m,ndim,np),taucld(m,ndim,np,2)
997
998 c-----output parameters
999
1000 _RL cc(m,n,3),tauclb(m,n,np),tauclf(m,n,np)
1001
1002 c-----temporary variables
1003
1004 integer i,j,k,im,it,ia,kk
1005 _RL fm,ft,fa,xai,taux
1006
1007 c-----pre-computed table
1008
1009 integer nm,nt,na
1010 parameter (nm=11,nt=9,na=11)
1011 _RL dm,dt,da,t1,caib(nm,nt,na),caif(nt,na)
1012 parameter (dm=0.1,dt=0.30103,da=0.1,t1=-0.9031)
1013
1014 c-----include the pre-computed table for cai
1015
1016 #include "cai-dat.h"
1017 save caib,caif
1018
1019
1020 c-----clouds within each of the high, middle, and low clouds are
1021 c assumed maximally overlapped, and the cloud cover (cc)
1022 c for a group is the maximum cloud cover of all the layers
1023 c in the group
1024
1025 do j=1,n
1026 do i=1,m
1027 cc(i,j,1)=0.0
1028 cc(i,j,2)=0.0
1029 cc(i,j,3)=0.0
1030 enddo
1031 enddo
1032
1033 do k=1,ict-1
1034 do j=1,n
1035 do i=1,m
1036 cc(i,j,1)=max(cc(i,j,1),fcld(i,j,k))
1037 enddo
1038 enddo
1039 enddo
1040
1041 do k=ict,icb-1
1042 do j=1,n
1043 do i=1,m
1044 cc(i,j,2)=max(cc(i,j,2),fcld(i,j,k))
1045 enddo
1046 enddo
1047 enddo
1048
1049 do k=icb,np
1050 do j=1,n
1051 do i=1,m
1052 cc(i,j,3)=max(cc(i,j,3),fcld(i,j,k))
1053 enddo
1054 enddo
1055 enddo
1056
1057 c-----scale the cloud optical thickness.
1058 c taucld(i,j,k,1) is the optical thickness for ice particles, and
1059 c taucld(i,j,k,2) is the optical thickness for liquid particles.
1060
1061 do k=1,np
1062
1063 if(k.lt.ict) then
1064 kk=1
1065 elseif(k.ge.ict .and. k.lt.icb) then
1066 kk=2
1067 else
1068 kk=3
1069 endif
1070
1071 do j=1,n
1072 do i=1,m
1073
1074 tauclb(i,j,k) = 0.0
1075 tauclf(i,j,k) = 0.0
1076
1077 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1078 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1079
1080 c-----normalize cloud cover
1081
1082 fa=fcld(i,j,k)/cc(i,j,kk)
1083
1084 c-----table look-up
1085
1086 taux=min(taux,32. _d 0)
1087
1088 fm=cosz(i,j)/dm
1089 ft=(log10(taux)-t1)/dt
1090 fa=fa/da
1091
1092 im=int(fm+1.5)
1093 it=int(ft+1.5)
1094 ia=int(fa+1.5)
1095
1096 im=max(im,2)
1097 it=max(it,2)
1098 ia=max(ia,2)
1099
1100 im=min(im,nm-1)
1101 it=min(it,nt-1)
1102 ia=min(ia,na-1)
1103
1104 fm=fm-float(im-1)
1105 ft=ft-float(it-1)
1106 fa=fa-float(ia-1)
1107
1108 c-----scale cloud optical thickness for beam radiation.
1109 c the scaling factor, xai, is a function of the solar zenith
1110 c angle, optical thickness, and cloud cover.
1111
1112 xai= (-caib(im-1,it,ia)*(1.-fm)+
1113 * caib(im+1,it,ia)*(1.+fm))*fm*.5+caib(im,it,ia)*(1.-fm*fm)
1114
1115 xai=xai+(-caib(im,it-1,ia)*(1.-ft)+
1116 * caib(im,it+1,ia)*(1.+ft))*ft*.5+caib(im,it,ia)*(1.-ft*ft)
1117
1118 xai=xai+(-caib(im,it,ia-1)*(1.-fa)+
1119 * caib(im,it,ia+1)*(1.+fa))*fa*.5+caib(im,it,ia)*(1.-fa*fa)
1120
1121 xai= xai-2.*caib(im,it,ia)
1122 xai=max(xai,0.0 _d 0)
1123
1124 tauclb(i,j,k) = taux*xai
1125
1126 c-----scale cloud optical thickness for diffuse radiation.
1127 c the scaling factor, xai, is a function of the cloud optical
1128 c thickness and cover but not the solar zenith angle.
1129
1130 xai= (-caif(it-1,ia)*(1.-ft)+
1131 * caif(it+1,ia)*(1.+ft))*ft*.5+caif(it,ia)*(1.-ft*ft)
1132
1133 xai=xai+(-caif(it,ia-1)*(1.-fa)+
1134 * caif(it,ia+1)*(1.+fa))*fa*.5+caif(it,ia)*(1.-fa*fa)
1135
1136 xai= xai-caif(it,ia)
1137 xai=max(xai,0.0 _d 0)
1138
1139 tauclf(i,j,k) = taux*xai
1140
1141 endif
1142
1143 enddo
1144 enddo
1145 enddo
1146
1147 return
1148 end
1149 c***********************************************************************
1150
1151 subroutine solir (m,n,ndim,np,wh,taucld,tauclb,tauclf,reff,
1152 * ict,icb,fcld,cc,taual,csm,rsirbm,rsirdf,flx,flc,fdirir,fdifir)
1153
1154 c************************************************************************
1155 c compute solar flux in the infrared region. The spectrum is divided
1156 c into three bands:
1157 c
1158 c band wavenumber(/cm) wavelength (micron)
1159 c 1 1000-4400 2.27-10.0
1160 c 2 4400-8200 1.22-2.27
1161 c 3 8200-14300 0.70-1.22
1162 c
1163 c----- Input parameters: units size
1164 c
1165 c number of soundings in zonal direction (m) n/d 1
1166 c number of soundings in meridional direction (n) n/d 1
1167 c maximum number of soundings in n/d 1
1168 c meridional direction (ndim)
1169 c number of atmospheric layers (np) n/d 1
1170 c layer water vapor content (wh) gm/cm^2 m*n*np
1171 c cloud optical thickness (taucld) n/d m*ndim*np*2
1172 c index 1 for ice paticles
1173 c index 2 for liquid particles
1174 c scaled cloud optical thickness n/d m*n*np
1175 c for beam radiation (tauclb)
1176 c scaled cloud optical thickness n/d m*n*np
1177 c for diffuse radiation (tauclf)
1178 c effective cloud-particle size (reff) micrometer m*ndim*np*2
1179 c index 1 for ice paticles
1180 c index 2 for liquid particles
1181 c level index separating high and n/d m*n
1182 c middle clouds (ict)
1183 c level index separating middle and n/d m*n
1184 c low clouds (icb)
1185 c cloud amount (fcld) fraction m*ndim*np
1186 c cloud amount of high, middle, and n/d m*n*3
1187 c low clouds (cc)
1188 c aerosol optical thickness (taual) n/d m*ndim*np
1189 c cosecant of the solar zenith angle (csm) n/d m*n
1190 c near ir surface albedo for beam fraction m*ndim
1191 c radiation (rsirbm)
1192 c near ir surface albedo for diffuse fraction m*ndim
1193 c radiation (rsirdf)
1194 c
1195 c----- output (updated) parameters:
1196 c
1197 c all-sky flux (downward-upward) (flx) fraction m*ndim*(np+1)
1198 c clear-sky flux (downward-upward) (flc) fraction m*ndim*(np+1)
1199 c all-sky direct downward ir flux at
1200 c the surface (fdirir) fraction m*ndim
1201 c all-sky diffuse downward ir flux at
1202 c the surface (fdifir) fraction m*ndim
1203 c
1204 c----- note: the following parameters must be specified by users:
1205 c aerosol single scattering albedo (ssaal) n/d nband
1206 c aerosol asymmetry factor (asyal) n/d nband
1207 c
1208 c*************************************************************************
1209
1210 implicit none
1211
1212 c-----Explicit Inline Directives
1213
1214 #ifdef CRAY
1215 #ifdef f77
1216 cfpp$ expand (deledd)
1217 cfpp$ expand (sagpol)
1218 cfpp$ expand (expmn)
1219 #endif
1220 #endif
1221 _RL expmn
1222
1223 c-----input parameters
1224
1225 integer m,n,ndim,np,ict,icb
1226 _RL taucld(m,ndim,np,2),reff(m,ndim,np,2),fcld(m,ndim,np)
1227 _RL tauclb(m,n,np),tauclf(m,n,np),cc(m,n,3)
1228 _RL rsirbm(m,ndim),rsirdf(m,ndim)
1229 _RL wh(m,n,np),taual(m,ndim,np),csm(m,n)
1230
1231 c-----output (updated) parameters
1232
1233 _RL flx(m,ndim,np+1),flc(m,ndim,np+1)
1234 _RL fdirir(m,ndim),fdifir(m,ndim)
1235
1236 c-----static parameters
1237
1238 integer nk,nband
1239 parameter (nk=10,nband=3)
1240 _RL xk(nk),hk(nband,nk),ssaal(nband),asyal(nband)
1241 _RL aia(nband,3),awa(nband,3),aig(nband,3),awg(nband,3)
1242
1243 c-----temporary array
1244
1245 integer ib,ik,i,j,k
1246 _RL ssacl(m,n,np),asycl(m,n,np)
1247 _RL rr(m,n,np+1,2),tt(m,n,np+1,2),td(m,n,np+1,2),
1248 * rs(m,n,np+1,2),ts(m,n,np+1,2)
1249 _RL fall(m,n,np+1),fclr(m,n,np+1)
1250 _RL fsdir(m,n),fsdif(m,n)
1251
1252 _RL tauwv,tausto,ssatau,asysto,tauto,ssato,asyto
1253 _RL taux,reff1,reff2,w1,w2,g1,g2
1254 _RL ssaclt(m,n),asyclt(m,n)
1255 _RL rr1t(m,n),tt1t(m,n),td1t(m,n),rs1t(m,n),ts1t(m,n)
1256 _RL rr2t(m,n),tt2t(m,n),td2t(m,n),rs2t(m,n),ts2t(m,n)
1257
1258 c-----water vapor absorption coefficient for 10 k-intervals.
1259 c unit: cm^2/gm
1260
1261 data xk/
1262 1 0.0010, 0.0133, 0.0422, 0.1334, 0.4217,
1263 2 1.334, 5.623, 31.62, 177.8, 1000.0/
1264
1265 c-----water vapor k-distribution function,
1266 c the sum of hk is 0.52926. unit: fraction
1267
1268 data hk/
1269 1 .01074,.08236,.20673, .00360,.01157,.03497,
1270 2 .00411,.01133,.03011, .00421,.01143,.02260,
1271 3 .00389,.01240,.01336, .00326,.01258,.00696,
1272 4 .00499,.01381,.00441, .00465,.00650,.00115,
1273 5 .00245,.00244,.00026, .00145,.00094,.00000/
1274
1275 c-----aerosol single-scattering albedo and asymmetry factor
1276
1277 data ssaal,asyal/nband*0.999,nband*0.850/
1278
1279 c-----coefficients for computing the single scattering albedo of
1280 c ice clouds from ssa=1-(aia(*,1)+aia(*,2)*reff+aia(*,3)*reff**2)
1281
1282 data aia/
1283 1 .08938331, .00215346,-.00000260,
1284 2 .00299387, .00073709, .00000746,
1285 3 -.00001038,-.00000134, .00000000/
1286
1287 c-----coefficients for computing the single scattering albedo of
1288 c liquid clouds from ssa=1-(awa(*,1)+awa(*,2)*reff+awa(*,3)*reff**2)
1289
1290 data awa/
1291 1 .01209318,-.00019934, .00000007,
1292 2 .01784739, .00088757, .00000845,
1293 3 -.00036910,-.00000650,-.00000004/
1294
1295 c-----coefficients for computing the asymmetry factor of ice clouds
1296 c from asycl=aig(*,1)+aig(*,2)*reff+aig(*,3)*reff**2
1297
1298 data aig/
1299 1 .84090400, .76098937, .74935228,
1300 2 .00126222, .00141864, .00119715,
1301 3 -.00000385,-.00000396,-.00000367/
1302
1303 c-----coefficients for computing the asymmetry factor of liquid clouds
1304 c from asycl=awg(*,1)+awg(*,2)*reff+awg(*,3)*reff**2
1305
1306 data awg/
1307 1 .83530748, .74513197, .79375035,
1308 2 .00257181, .01370071, .00832441,
1309 3 .00005519,-.00038203,-.00023263/
1310
1311 c-----initialize surface fluxes, reflectances, and transmittances
1312
1313 do j= 1, n
1314 do i= 1, m
1315 fdirir(i,j)=0.0
1316 fdifir(i,j)=0.0
1317 rr(i,j,np+1,1)=rsirbm(i,j)
1318 rr(i,j,np+1,2)=rsirbm(i,j)
1319 rs(i,j,np+1,1)=rsirdf(i,j)
1320 rs(i,j,np+1,2)=rsirdf(i,j)
1321 td(i,j,np+1,1)=0.0
1322 td(i,j,np+1,2)=0.0
1323 tt(i,j,np+1,1)=0.0
1324 tt(i,j,np+1,2)=0.0
1325 ts(i,j,np+1,1)=0.0
1326 ts(i,j,np+1,2)=0.0
1327 enddo
1328 enddo
1329
1330 c-----integration over spectral bands
1331
1332 do 100 ib=1,nband
1333
1334 c-----compute cloud single scattering albedo and asymmetry factor
1335 c for a mixture of ice and liquid particles.
1336
1337 do k= 1, np
1338
1339 do j= 1, n
1340 do i= 1, m
1341
1342 ssaclt(i,j)=1.0
1343 asyclt(i,j)=1.0
1344
1345 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1346 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1347
1348 reff1=min(reff(i,j,k,1),130. _d 0)
1349 reff2=min(reff(i,j,k,2),20.0 _d 0)
1350
1351 w1=(1.-(aia(ib,1)+(aia(ib,2)+
1352 * aia(ib,3)*reff1)*reff1))*taucld(i,j,k,1)
1353 w2=(1.-(awa(ib,1)+(awa(ib,2)+
1354 * awa(ib,3)*reff2)*reff2))*taucld(i,j,k,2)
1355 ssaclt(i,j)=(w1+w2)/taux
1356
1357 g1=(aig(ib,1)+(aig(ib,2)+aig(ib,3)*reff1)*reff1)*w1
1358 g2=(awg(ib,1)+(awg(ib,2)+awg(ib,3)*reff2)*reff2)*w2
1359 asyclt(i,j)=(g1+g2)/(w1+w2)
1360
1361 endif
1362
1363 enddo
1364 enddo
1365
1366 do j=1,n
1367 do i=1,m
1368 ssacl(i,j,k)=ssaclt(i,j)
1369 enddo
1370 enddo
1371 do j=1,n
1372 do i=1,m
1373 asycl(i,j,k)=asyclt(i,j)
1374 enddo
1375 enddo
1376
1377 enddo
1378
1379 c-----integration over the k-distribution function
1380
1381 do 200 ik=1,nk
1382
1383 do 300 k= 1, np
1384
1385 do j= 1, n
1386 do i= 1, m
1387
1388 tauwv=xk(ik)*wh(i,j,k)
1389
1390 c-----compute total optical thickness, single scattering albedo,
1391 c and asymmetry factor.
1392
1393 tausto=tauwv+taual(i,j,k)+1.0e-8
1394 ssatau=ssaal(ib)*taual(i,j,k)
1395 asysto=asyal(ib)*ssaal(ib)*taual(i,j,k)
1396
1397 c-----compute reflectance and transmittance for cloudless layers
1398
1399 tauto=tausto
1400 ssato=ssatau/tauto+1.0e-8
1401
1402 c if (ssato .gt. 0.001) then
1403
1404 c ssato=min(ssato,0.999999 _d 0)
1405 c asyto=asysto/(ssato*tauto)
1406
1407 c call deledd(tauto,ssato,asyto,csm(i,j),
1408 c * rr1t(i,j),tt1t(i,j),td1t(i,j))
1409
1410 c call sagpol (tauto,ssato,asyto,rs1t(i,j),ts1t(i,j))
1411
1412 c else
1413
1414 td1t(i,j)=expmn(-tauto*csm(i,j))
1415 ts1t(i,j)=expmn(-1.66*tauto)
1416 tt1t(i,j)=0.0
1417 rr1t(i,j)=0.0
1418 rs1t(i,j)=0.0
1419
1420 c endif
1421
1422 c-----compute reflectance and transmittance for cloud layers
1423
1424 if (tauclb(i,j,k) .lt. 0.01) then
1425
1426 rr2t(i,j)=rr1t(i,j)
1427 tt2t(i,j)=tt1t(i,j)
1428 td2t(i,j)=td1t(i,j)
1429 rs2t(i,j)=rs1t(i,j)
1430 ts2t(i,j)=ts1t(i,j)
1431
1432 else
1433
1434 tauto=tausto+tauclb(i,j,k)
1435 ssato=(ssatau+ssacl(i,j,k)*tauclb(i,j,k))/tauto+1.0e-8
1436 ssato=min(ssato,0.999999 _d 0)
1437 asyto=(asysto+asycl(i,j,k)*ssacl(i,j,k)*tauclb(i,j,k))/
1438 * (ssato*tauto)
1439
1440 call deledd(tauto,ssato,asyto,csm(i,j),
1441 * rr2t(i,j),tt2t(i,j),td2t(i,j))
1442
1443 tauto=tausto+tauclf(i,j,k)
1444 ssato=(ssatau+ssacl(i,j,k)*tauclf(i,j,k))/tauto+1.0e-8
1445 ssato=min(ssato,0.999999 _d 0)
1446 asyto=(asysto+asycl(i,j,k)*ssacl(i,j,k)*tauclf(i,j,k))/
1447 * (ssato*tauto)
1448
1449 call sagpol (tauto,ssato,asyto,rs2t(i,j),ts2t(i,j))
1450
1451 endif
1452
1453 enddo
1454 enddo
1455
1456 c-----the following code appears to be awkward, but it is efficient
1457 c in certain parallel processors.
1458
1459 do j=1,n
1460 do i=1,m
1461 rr(i,j,k,1)=rr1t(i,j)
1462 enddo
1463 enddo
1464 do j=1,n
1465 do i=1,m
1466 tt(i,j,k,1)=tt1t(i,j)
1467 enddo
1468 enddo
1469 do j=1,n
1470 do i=1,m
1471 td(i,j,k,1)=td1t(i,j)
1472 enddo
1473 enddo
1474 do j=1,n
1475 do i=1,m
1476 rs(i,j,k,1)=rs1t(i,j)
1477 enddo
1478 enddo
1479 do j=1,n
1480 do i=1,m
1481 ts(i,j,k,1)=ts1t(i,j)
1482 enddo
1483 enddo
1484
1485 do j=1,n
1486 do i=1,m
1487 rr(i,j,k,2)=rr2t(i,j)
1488 enddo
1489 enddo
1490 do j=1,n
1491 do i=1,m
1492 tt(i,j,k,2)=tt2t(i,j)
1493 enddo
1494 enddo
1495 do j=1,n
1496 do i=1,m
1497 td(i,j,k,2)=td2t(i,j)
1498 enddo
1499 enddo
1500 do j=1,n
1501 do i=1,m
1502 rs(i,j,k,2)=rs2t(i,j)
1503 enddo
1504 enddo
1505 do j=1,n
1506 do i=1,m
1507 ts(i,j,k,2)=ts2t(i,j)
1508 enddo
1509 enddo
1510
1511 300 continue
1512
1513 c-----flux calculations
1514
1515 do k= 1, np+1
1516 do j= 1, n
1517 do i= 1, m
1518 fclr(i,j,k) = 0.
1519 fall(i,j,k) = 0.
1520 enddo
1521 enddo
1522 enddo
1523 do j= 1, n
1524 do i= 1, m
1525 fsdir(i,j) = 0.
1526 fsdif(i,j) = 0.
1527 enddo
1528 enddo
1529
1530 call cldflx (m,n,np,ict,icb,cc,rr,tt,td,rs,ts,
1531 * fclr,fall,fsdir,fsdif)
1532
1533 do k= 1, np+1
1534 do j= 1, n
1535 do i= 1, m
1536 flx(i,j,k) = flx(i,j,k)+fall(i,j,k)*hk(ib,ik)
1537 enddo
1538 enddo
1539 do j= 1, n
1540 do i= 1, m
1541 flc(i,j,k) = flc(i,j,k)+fclr(i,j,k)*hk(ib,ik)
1542 enddo
1543 enddo
1544 enddo
1545
1546 do j= 1, n
1547 do i= 1, m
1548 fdirir(i,j) = fdirir(i,j)+fsdir(i,j)*hk(ib,ik)
1549 fdifir(i,j) = fdifir(i,j)+fsdif(i,j)*hk(ib,ik)
1550 enddo
1551 enddo
1552
1553 200 continue
1554 100 continue
1555
1556 return
1557 end
1558
1559 c************************************************************************
1560
1561 subroutine soluv (m,n,ndim,np,oh,dp,taucld,tauclb,tauclf,reff,
1562 * ict,icb,fcld,cc,taual,csm,rsuvbm,rsuvdf,flx,flc
1563 * ,fdirpar,fdifpar,fdiruv,fdifuv)
1564
1565 c************************************************************************
1566 c compute solar fluxes in the uv+visible region. the spectrum is
1567 c grouped into 8 bands:
1568 c
1569 c Band Micrometer
1570 c
1571 c UV-C 1. .175 - .225
1572 c 2. .225 - .245
1573 c .260 - .280
1574 c 3. .245 - .260
1575 c
1576 c UV-B 4. .280 - .295
1577 c 5. .295 - .310
1578 c 6. .310 - .320
1579 c
1580 c UV-A 7. .320 - .400
1581 c
1582 c PAR 8. .400 - .700
1583 c
1584 c----- Input parameters: units size
1585 c
1586 c number of soundings in zonal direction (m) n/d 1
1587 c number of soundings in meridional direction (n) n/d 1
1588 c maximum number of soundings in n/d 1
1589 c meridional direction (ndim)
1590 c number of atmospheric layers (np) n/d 1
1591 c layer ozone content (oh) (cm-atm)stp m*n*np
1592 c layer pressure thickness (dp) mb m*n*np
1593 c cloud optical thickness (taucld) n/d m*ndim*np*2
1594 c index 1 for ice paticles
1595 c index 2 for liquid particles
1596 c scaled cloud optical thickness n/d m*n*np
1597 c for beam radiation (tauclb)
1598 c scaled cloud optical thickness n/d m*n*np
1599 c for diffuse radiation (tauclf)
1600 c effective cloud-particle size (reff) micrometer m*ndim*np*2
1601 c index 1 for ice paticles
1602 c index 2 for liquid particles
1603 c level indiex separating high and n/d m*n
1604 c middle clouds (ict)
1605 c level indiex separating middle and n/d m*n
1606 c low clouds (icb)
1607 c cloud amount (fcld) fraction m*ndim*np
1608 c cloud amounts of high, middle, and n/d m*n*3
1609 c low clouds (cc)
1610 c aerosol optical thickness (taual) n/d m*ndim*np
1611 c cosecant of the solar zenith angle (csm) n/d m*n
1612 c uv+par surface albedo for beam fraction m*ndim
1613 c radiation (rsuvbm)
1614 c uv+par surface albedo for diffuse fraction m*ndim
1615 c radiation (rsuvdf)
1616 c
1617 c----- output (updated) parameters:
1618 c
1619 c all-sky net downward flux (flx) fraction m*ndim*(np+1)
1620 c clear-sky net downward flux (flc) fraction m*ndim*(np+1)
1621 c all-sky direct downward par flux at
1622 c the surface (fdirpar) fraction m*ndim
1623 c all-sky diffuse downward par flux at
1624 c the surface (fdifpar) fraction m*ndim
1625 c
1626 c----- note: the following parameters must be specified by users:
1627 c
1628 c aerosol single scattering albedo (ssaal) n/d 1
1629 c aerosol asymmetry factor (asyal) n/d 1
1630 c
1631 *
1632 c***********************************************************************
1633
1634 implicit none
1635
1636 c-----Explicit Inline Directives
1637
1638 #ifdef CRAY
1639 #ifdef f77
1640 cfpp$ expand (deledd)
1641 cfpp$ expand (sagpol)
1642 #endif
1643 #endif
1644
1645 c-----input parameters
1646
1647 integer m,n,ndim,np,ict,icb
1648 _RL taucld(m,ndim,np,2),reff(m,ndim,np,2),fcld(m,ndim,np)
1649 _RL tauclb(m,n,np),tauclf(m,n,np),cc(m,n,3)
1650 _RL oh(m,n,np),dp(m,n,np),taual(m,ndim,np)
1651 _RL rsuvbm(m,ndim),rsuvdf(m,ndim),csm(m,n)
1652
1653 c-----output (updated) parameter
1654
1655 _RL flx(m,ndim,np+1),flc(m,ndim,np+1)
1656 _RL fdirpar(m,ndim),fdifpar(m,ndim)
1657 _RL fdiruv(m,ndim),fdifuv(m,ndim)
1658
1659 c-----static parameters
1660
1661 integer nband
1662 parameter (nband=8)
1663 _RL hk(nband),xk(nband),ry(nband)
1664 _RL asyal(nband),ssaal(nband),aig(3),awg(3)
1665
1666 c-----temporary array
1667
1668 integer i,j,k,ib
1669 _RL taurs,tauoz,tausto,ssatau,asysto,tauto,ssato,asyto
1670 _RL taux,reff1,reff2,g1,g2,asycl(m,n,np)
1671 _RL td(m,n,np+1,2),rr(m,n,np+1,2),tt(m,n,np+1,2),
1672 * rs(m,n,np+1,2),ts(m,n,np+1,2)
1673 _RL fall(m,n,np+1),fclr(m,n,np+1),fsdir(m,n),fsdif(m,n)
1674 _RL asyclt(m,n)
1675 _RL rr1t(m,n),tt1t(m,n),td1t(m,n),rs1t(m,n),ts1t(m,n)
1676 _RL rr2t(m,n),tt2t(m,n),td2t(m,n),rs2t(m,n),ts2t(m,n)
1677
1678 c-----hk is the fractional extra-terrestrial solar flux.
1679 c the sum of hk is 0.47074.
1680
1681 data hk/.00057, .00367, .00083, .00417,
1682 * .00600, .00556, .05913, .39081/
1683
1684 c-----xk is the ozone absorption coefficient. unit: /(cm-atm)stp
1685
1686 data xk /30.47, 187.2, 301.9, 42.83,
1687 * 7.09, 1.25, 0.0345, 0.0539/
1688
1689 c-----ry is the extinction coefficient for Rayleigh scattering.
1690 c unit: /mb.
1691
1692 data ry /.00604, .00170, .00222, .00132,
1693 * .00107, .00091, .00055, .00012/
1694
1695 c-----aerosol single-scattering albedo and asymmetry factor
1696
1697 data ssaal,asyal/nband*0.999,nband*0.850/
1698
1699 c-----coefficients for computing the asymmetry factor of ice clouds
1700 c from asycl=aig(*,1)+aig(*,2)*reff+aig(*,3)*reff**2
1701
1702 data aig/.74625000,.00105410,-.00000264/
1703
1704 c-----coefficients for computing the asymmetry factor of liquid
1705 c clouds from asycl=awg(*,1)+awg(*,2)*reff+awg(*,3)*reff**2
1706
1707 data awg/.82562000,.00529000,-.00014866/
1708
1709 c-----initialize surface reflectances and transmittances
1710
1711 do j= 1, n
1712 do i= 1, m
1713 rr(i,j,np+1,1)=rsuvbm(i,j)
1714 rr(i,j,np+1,2)=rsuvbm(i,j)
1715 rs(i,j,np+1,1)=rsuvdf(i,j)
1716 rs(i,j,np+1,2)=rsuvdf(i,j)
1717 td(i,j,np+1,1)=0.0
1718 td(i,j,np+1,2)=0.0
1719 tt(i,j,np+1,1)=0.0
1720 tt(i,j,np+1,2)=0.0
1721 ts(i,j,np+1,1)=0.0
1722 ts(i,j,np+1,2)=0.0
1723 enddo
1724 enddo
1725
1726 c-----compute cloud asymmetry factor for a mixture of
1727 c liquid and ice particles. unit of reff is micrometers.
1728
1729 do k= 1, np
1730
1731 do j= 1, n
1732 do i= 1, m
1733
1734 asyclt(i,j)=1.0
1735
1736 taux=taucld(i,j,k,1)+taucld(i,j,k,2)
1737 if (taux.gt.0.05 .and. fcld(i,j,k).gt.0.01) then
1738
1739 reff1=min(reff(i,j,k,1),130. _d 0)
1740 reff2=min(reff(i,j,k,2),20.0 _d 0)
1741
1742 g1=(aig(1)+(aig(2)+aig(3)*reff1)*reff1)*taucld(i,j,k,1)
1743 g2=(awg(1)+(awg(2)+awg(3)*reff2)*reff2)*taucld(i,j,k,2)
1744 asyclt(i,j)=(g1+g2)/taux
1745
1746 endif
1747
1748 enddo
1749 enddo
1750
1751 do j=1,n
1752 do i=1,m
1753 asycl(i,j,k)=asyclt(i,j)
1754 enddo
1755 enddo
1756
1757 enddo
1758
1759 do j=1,n
1760 do i=1,m
1761 fdiruv(i,j)=0.0
1762 fdifuv(i,j)=0.0
1763 enddo
1764 enddo
1765
1766 c-----integration over spectral bands
1767
1768 do 100 ib=1,nband
1769
1770 do 300 k= 1, np
1771
1772 do j= 1, n
1773 do i= 1, m
1774
1775 c-----compute ozone and rayleigh optical thicknesses
1776
1777 taurs=ry(ib)*dp(i,j,k)
1778 tauoz=xk(ib)*oh(i,j,k)
1779
1780 c-----compute total optical thickness, single scattering albedo,
1781 c and asymmetry factor
1782
1783 tausto=taurs+tauoz+taual(i,j,k)+1.0e-8
1784 ssatau=ssaal(ib)*taual(i,j,k)+taurs
1785 asysto=asyal(ib)*ssaal(ib)*taual(i,j,k)
1786
1787 c-----compute reflectance and transmittance for cloudless layers
1788
1789 tauto=tausto
1790 ssato=ssatau/tauto+1.0e-8
1791 ssato=min(ssato,0.999999 _d 0)
1792 asyto=asysto/(ssato*tauto)
1793
1794 call deledd(tauto,ssato,asyto,csm(i,j),
1795 * rr1t(i,j),tt1t(i,j),td1t(i,j))
1796
1797 call sagpol (tauto,ssato,asyto,rs1t(i,j),ts1t(i,j))
1798
1799 c-----compute reflectance and transmittance for cloud layers
1800
1801 if (tauclb(i,j,k) .lt. 0.01) then
1802
1803 rr2t(i,j)=rr1t(i,j)
1804 tt2t(i,j)=tt1t(i,j)
1805 td2t(i,j)=td1t(i,j)
1806 rs2t(i,j)=rs1t(i,j)
1807 ts2t(i,j)=ts1t(i,j)
1808
1809 else
1810
1811 tauto=tausto+tauclb(i,j,k)
1812 ssato=(ssatau+tauclb(i,j,k))/tauto+1.0e-8
1813 ssato=min(ssato,0.999999 _d 0)
1814 asyto=(asysto+asycl(i,j,k)*tauclb(i,j,k))/(ssato*tauto)
1815
1816 call deledd(tauto,ssato,asyto,csm(i,j),
1817 * rr2t(i,j),tt2t(i,j),td2t(i,j))
1818
1819 tauto=tausto+tauclf(i,j,k)
1820 ssato=(ssatau+tauclf(i,j,k))/tauto+1.0e-8
1821 ssato=min(ssato,0.999999 _d 0)
1822 asyto=(asysto+asycl(i,j,k)*tauclf(i,j,k))/(ssato*tauto)
1823
1824 call sagpol (tauto,ssato,asyto,rs2t(i,j),ts2t(i,j))
1825
1826 endif
1827
1828 enddo
1829 enddo
1830
1831 do j=1,n
1832 do i=1,m
1833 rr(i,j,k,1)=rr1t(i,j)
1834 enddo
1835 enddo
1836 do j=1,n
1837 do i=1,m
1838 tt(i,j,k,1)=tt1t(i,j)
1839 enddo
1840 enddo
1841 do j=1,n
1842 do i=1,m
1843 td(i,j,k,1)=td1t(i,j)
1844 enddo
1845 enddo
1846 do j=1,n
1847 do i=1,m
1848 rs(i,j,k,1)=rs1t(i,j)
1849 enddo
1850 enddo
1851 do j=1,n
1852 do i=1,m
1853 ts(i,j,k,1)=ts1t(i,j)
1854 enddo
1855 enddo
1856
1857 do j=1,n
1858 do i=1,m
1859 rr(i,j,k,2)=rr2t(i,j)
1860 enddo
1861 enddo
1862 do j=1,n
1863 do i=1,m
1864 tt(i,j,k,2)=tt2t(i,j)
1865 enddo
1866 enddo
1867 do j=1,n
1868 do i=1,m
1869 td(i,j,k,2)=td2t(i,j)
1870 enddo
1871 enddo
1872 do j=1,n
1873 do i=1,m
1874 rs(i,j,k,2)=rs2t(i,j)
1875 enddo
1876 enddo
1877 do j=1,n
1878 do i=1,m
1879 ts(i,j,k,2)=ts2t(i,j)
1880 enddo
1881 enddo
1882
1883 300 continue
1884
1885 c-----flux calculations
1886
1887 do k= 1, np+1
1888 do j= 1, n
1889 do i= 1, m
1890 fclr(i,j,k) = 0.
1891 fall(i,j,k) = 0.
1892 enddo
1893 enddo
1894 enddo
1895 do j= 1, n
1896 do i= 1, m
1897 fsdir(i,j) = 0.
1898 fsdif(i,j) = 0.
1899 enddo
1900 enddo
1901 call cldflx (m,n,np,ict,icb,cc,rr,tt,td,rs,ts,
1902 * fclr,fall,fsdir,fsdif)
1903
1904 do k= 1, np+1
1905 do j= 1, n
1906 do i= 1, m
1907 flx(i,j,k)=flx(i,j,k)+fall(i,j,k)*hk(ib)
1908 enddo
1909 enddo
1910 do j= 1, n
1911 do i= 1, m
1912 flc(i,j,k)=flc(i,j,k)+fclr(i,j,k)*hk(ib)
1913 enddo
1914 enddo
1915 enddo
1916
1917 if(ib.eq.nband) then
1918 do j=1,n
1919 do i=1,m
1920 fdirpar(i,j)=fsdir(i,j)*hk(ib)
1921 fdifpar(i,j)=fsdif(i,j)*hk(ib)
1922 enddo
1923 enddo
1924 else
1925 do j=1,n
1926 do i=1,m
1927 fdiruv(i,j)=fdiruv(i,j)+fsdir(i,j)*hk(ib)
1928 fdifuv(i,j)=fdifuv(i,j)+fsdif(i,j)*hk(ib)
1929 enddo
1930 enddo
1931 endif
1932
1933 100 continue
1934
1935 return
1936 end
1937
1938 c*********************************************************************
1939
1940 subroutine deledd(tau,ssc,g0,csm,rr,tt,td)
1941
1942 c*********************************************************************
1943 c
1944 c-----uses the delta-eddington approximation to compute the
1945 c bulk scattering properties of a single layer
1946 c coded following King and Harshvardhan (JAS, 1986)
1947 c
1948 c inputs:
1949 c
1950 c tau: the effective optical thickness
1951 c ssc: the effective single scattering albedo
1952 c g0: the effective asymmetry factor
1953 c csm: the effective secant of the zenith angle
1954 c
1955 c outputs:
1956 c
1957 c rr: the layer reflection of the direct beam
1958 c tt: the layer diffuse transmission of the direct beam
1959 c td: the layer direct transmission of the direct beam
1960
1961 c*********************************************************************
1962
1963 implicit none
1964
1965 c-----Explicit Inline Directives
1966
1967 #ifdef CRAY
1968 #ifdef f77
1969 cfpp$ expand (expmn)
1970 #endif
1971 #endif
1972 _RL expmn
1973
1974 _RL zero,one,two,three,four,fourth,seven,tumin
1975 parameter (one=1., three=3.)
1976 parameter (seven=7., two=2.)
1977 parameter (four=4., fourth=.25)
1978 parameter (zero=0., tumin=1.e-20)
1979
1980 c-----input parameters
1981 _RL tau,ssc,g0,csm
1982
1983 c-----output parameters
1984 _RL rr,tt,td
1985
1986 c-----temporary parameters
1987
1988 _RL zth,ff,xx,taup,sscp,gp,gm1,gm2,gm3,akk,alf1,alf2,
1989 * all1,bll,st7,st8,cll,dll,fll,ell,st1,st2,st3,st4
1990 c
1991 zth = one / csm
1992
1993 c delta-eddington scaling of single scattering albedo,
1994 c optical thickness, and asymmetry factor,
1995 c K & H eqs(27-29)
1996
1997 ff = g0*g0
1998 xx = one-ff*ssc
1999 taup= tau*xx
2000 sscp= ssc*(one-ff)/xx
2001 gp = g0/(one+g0)
2002
2003 c extinction of the direct beam transmission
2004
2005 td = expmn(-taup*csm)
2006
2007 c gamma1, gamma2, gamma3 and gamma4. see table 2 and eq(26) K & H
2008 c ssc and gp are the d-s single scattering
2009 c albedo and asymmetry factor.
2010
2011 xx = three*gp
2012 gm1 = (seven - sscp*(four+xx))*fourth
2013 gm2 = -(one - sscp*(four-xx))*fourth
2014 gm3 = (two - zth*xx )*fourth
2015
2016 c akk is k as defined in eq(25) of K & H
2017
2018 xx = gm1-gm2
2019 akk = sqrt((gm1+gm2)*xx)
2020
2021 c alf1 and alf2 are alpha1 and alpha2 from eqs (23) & (24) of K & H
2022
2023 alf1 = gm1 - gm3 * xx
2024 alf2 = gm2 + gm3 * xx
2025
2026 c all1 is last term in eq(21) of K & H
2027 c bll is last term in eq(22) of K & H
2028
2029 xx = akk * two
2030 all1 = (gm3 - alf2 * zth )*xx*td
2031 bll = (one - gm3 + alf1*zth)*xx
2032
2033 xx = akk * zth
2034 st7 = one - xx
2035 st8 = one + xx
2036
2037 xx = akk * gm3
2038 cll = (alf2 + xx) * st7
2039 dll = (alf2 - xx) * st8
2040
2041 xx = akk * (one-gm3)
2042 fll = (alf1 + xx) * st8
2043 ell = (alf1 - xx) * st7
2044
2045 st3 = max(abs(st7*st8),tumin)
2046 st3 = sign (st3,st7)
2047
2048 st2 = expmn(-akk*taup)
2049 st4 = st2 * st2
2050
2051 st1 = sscp / ((akk+gm1 + (akk-gm1)*st4) * st3)
2052
2053 c rr is r-hat of eq(21) of K & H
2054 c tt is diffuse part of t-hat of eq(22) of K & H
2055
2056 rr = ( cll-dll*st4 -all1*st2)*st1
2057 tt = - ((fll-ell*st4)*td-bll*st2)*st1
2058
2059 rr = max(rr,zero)
2060 tt = max(tt,zero)
2061
2062 return
2063 end
2064
2065 c*********************************************************************
2066
2067 subroutine sagpol(tau,ssc,g0,rll,tll)
2068
2069 c*********************************************************************
2070 c-----transmittance (tll) and reflectance (rll) of diffuse radiation
2071 c follows Sagan and Pollock (JGR, 1967).
2072 c also, eq.(31) of Lacis and Hansen (JAS, 1974).
2073 c
2074 c-----input parameters:
2075 c
2076 c tau: the effective optical thickness
2077 c ssc: the effective single scattering albedo
2078 c g0: the effective asymmetry factor
2079 c
2080 c-----output parameters:
2081 c
2082 c rll: the layer reflection of diffuse radiation
2083 c tll: the layer transmission of diffuse radiation
2084 c
2085 c*********************************************************************
2086
2087 implicit none
2088
2089 c-----Explicit Inline Directives
2090
2091 #ifdef CRAY
2092 #ifdef f77
2093 cfpp$ expand (expmn)
2094 #endif
2095 #endif
2096 _RL expmn
2097
2098 _RL one,three,four
2099 parameter (one=1., three=3., four=4.)
2100
2101 c-----output parameters:
2102
2103 _RL tau,ssc,g0
2104
2105 c-----output parameters:
2106
2107 _RL rll,tll
2108
2109 c-----temporary arrays
2110
2111 _RL xx,uuu,ttt,emt,up1,um1,st1
2112
2113 xx = one-ssc*g0
2114 uuu = sqrt( xx/(one-ssc))
2115 ttt = sqrt( xx*(one-ssc)*three )*tau
2116 emt = expmn(-ttt)
2117 up1 = uuu + one
2118 um1 = uuu - one
2119 xx = um1*emt
2120 st1 = one / ((up1+xx) * (up1-xx))
2121 rll = up1*um1*(one-emt*emt)*st1
2122 tll = uuu*four*emt *st1
2123
2124 return
2125 end
2126
2127 c*******************************************************************
2128
2129 function expmn(fin)
2130
2131 c*******************************************************************
2132 c compute exponential for arguments in the range 0> fin > -10.
2133 c*******************************************************************
2134 implicit none
2135 _RL fin,expmn
2136
2137 _RL one,expmin,e1,e2,e3,e4
2138 parameter (one=1.0, expmin=-10.0)
2139 parameter (e1=1.0, e2=-2.507213e-1)
2140 parameter (e3=2.92732e-2, e4=-3.827800e-3)
2141
2142 if (fin .lt. expmin) fin = expmin
2143 expmn = ((e4*fin + e3)*fin+e2)*fin+e1
2144 expmn = expmn * expmn
2145 expmn = one / (expmn * expmn)
2146
2147 return
2148 end
2149
2150 c*******************************************************************
2151
2152 subroutine cldflx (m,n,np,ict,icb,cc,rr,tt,td,rs,ts,
2153 * fclr,fall,fsdir,fsdif)
2154
2155 c*******************************************************************
2156 c compute upward and downward fluxes using a two-stream adding method
2157 c following equations (3)-(5) of Chou (1992, JAS).
2158 c
2159 c clouds are grouped into high, middle, and low clouds which are
2160 c assumed randomly overlapped. It involves eight sets of calculations.
2161 c In each set of calculations, each atmospheric layer is homogeneous,
2162 c either with clouds or without clouds.
2163
2164 c input parameters:
2165 c m: number of soundings in zonal direction
2166 c n: number of soundings in meridional direction
2167 c np: number of atmospheric layers
2168 c ict: the level separating high and middle clouds
2169 c icb: the level separating middle and low clouds
2170 c cc: effective cloud covers for high, middle and low clouds
2171 c tt: diffuse transmission of a layer illuminated by beam radiation
2172 c td: direct beam tranmssion
2173 c ts: transmission of a layer illuminated by diffuse radiation
2174 c rr: reflection of a layer illuminated by beam radiation
2175 c rs: reflection of a layer illuminated by diffuse radiation
2176 c
2177 c output parameters:
2178 c fclr: clear-sky flux (downward minus upward)
2179 c fall: all-sky flux (downward minus upward)
2180 c fdndir: surface direct downward flux
2181 c fdndif: surface diffuse downward flux
2182 c*********************************************************************c
2183
2184 implicit none
2185
2186 c-----input parameters
2187
2188 integer m,n,np,ict,icb
2189
2190 _RL rr(m,n,np+1,2),tt(m,n,np+1,2),td(m,n,np+1,2)
2191 _RL rs(m,n,np+1,2),ts(m,n,np+1,2)
2192 _RL cc(m,n,3)
2193
2194 c-----temporary array
2195
2196 integer i,j,k,ih,im,is
2197 _RL rra(m,n,np+1,2,2),tta(m,n,np+1,2,2),tda(m,n,np+1,2,2)
2198 _RL rsa(m,n,np+1,2,2),rxa(m,n,np+1,2,2)
2199 _RL ch(m,n),cm(m,n),ct(m,n),flxdn(m,n,np+1)
2200 _RL fdndir(m,n),fdndif(m,n),fupdif
2201 _RL denm,xx
2202
2203 c-----output parameters
2204
2205 _RL fclr(m,n,np+1),fall(m,n,np+1)
2206 _RL fsdir(m,n),fsdif(m,n)
2207
2208 c-----initialize all-sky flux (fall) and surface downward fluxes
2209
2210 do k=1,np+1
2211 do j=1,n
2212 do i=1,m
2213 fall(i,j,k)=0.0
2214 enddo
2215 enddo
2216 enddo
2217
2218 do j=1,n
2219 do i=1,m
2220 fsdir(i,j)=0.0
2221 fsdif(i,j)=0.0
2222 enddo
2223 enddo
2224
2225 c-----compute transmittances and reflectances for a composite of
2226 c layers. layers are added one at a time, going down from the top.
2227 c tda is the composite transmittance illuminated by beam radiation
2228 c tta is the composite diffuse transmittance illuminated by
2229 c beam radiation
2230 c rsa is the composite reflectance illuminated from below
2231 c by diffuse radiation
2232 c tta and rsa are computed from eqs. (4b) and (3b) of Chou
2233
2234 c-----for high clouds. indices 1 and 2 denote clear and cloudy
2235 c situations, respectively.
2236
2237 do 10 ih=1,2
2238
2239 do j= 1, n
2240 do i= 1, m
2241 tda(i,j,1,ih,1)=td(i,j,1,ih)
2242 tta(i,j,1,ih,1)=tt(i,j,1,ih)
2243 rsa(i,j,1,ih,1)=rs(i,j,1,ih)
2244 tda(i,j,1,ih,2)=td(i,j,1,ih)
2245 tta(i,j,1,ih,2)=tt(i,j,1,ih)
2246 rsa(i,j,1,ih,2)=rs(i,j,1,ih)
2247 enddo
2248 enddo
2249
2250 do k= 2, ict-1
2251 do j= 1, n
2252 do i= 1, m
2253 denm = ts(i,j,k,ih)/( 1.-rsa(i,j,k-1,ih,1)*rs(i,j,k,ih))
2254 tda(i,j,k,ih,1)= tda(i,j,k-1,ih,1)*td(i,j,k,ih)
2255 tta(i,j,k,ih,1)= tda(i,j,k-1,ih,1)*tt(i,j,k,ih)
2256 * +(tda(i,j,k-1,ih,1)*rr(i,j,k,ih)
2257 * *rsa(i,j,k-1,ih,1)+tta(i,j,k-1,ih,1))*denm
2258 rsa(i,j,k,ih,1)= rs(i,j,k,ih)+ts(i,j,k,ih)
2259 * *rsa(i,j,k-1,ih,1)*denm
2260 tda(i,j,k,ih,2)= tda(i,j,k,ih,1)
2261 tta(i,j,k,ih,2)= tta(i,j,k,ih,1)
2262 rsa(i,j,k,ih,2)= rsa(i,j,k,ih,1)
2263 enddo
2264 enddo
2265 enddo
2266
2267 c-----for middle clouds
2268
2269 do 10 im=1,2
2270
2271 do k= ict, icb-1
2272 do j= 1, n
2273 do i= 1, m
2274 denm = ts(i,j,k,im)/( 1.-rsa(i,j,k-1,ih,im)*rs(i,j,k,im))
2275 tda(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*td(i,j,k,im)
2276 tta(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*tt(i,j,k,im)
2277 * +(tda(i,j,k-1,ih,im)*rr(i,j,k,im)
2278 * *rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2279 rsa(i,j,k,ih,im)= rs(i,j,k,im)+ts(i,j,k,im)
2280 * *rsa(i,j,k-1,ih,im)*denm
2281 enddo
2282 enddo
2283 enddo
2284
2285 10 continue
2286
2287 c-----layers are added one at a time, going up from the surface.
2288 c rra is the composite reflectance illuminated by beam radiation
2289 c rxa is the composite reflectance illuminated from above
2290 c by diffuse radiation
2291 c rra and rxa are computed from eqs. (4a) and (3a) of Chou
2292
2293 c-----for the low clouds
2294
2295 do 20 is=1,2
2296
2297 do j= 1, n
2298 do i= 1, m
2299 rra(i,j,np+1,1,is)=rr(i,j,np+1,is)
2300 rxa(i,j,np+1,1,is)=rs(i,j,np+1,is)
2301 rra(i,j,np+1,2,is)=rr(i,j,np+1,is)
2302 rxa(i,j,np+1,2,is)=rs(i,j,np+1,is)
2303 enddo
2304 enddo
2305
2306 do k=np,icb,-1
2307 do j= 1, n
2308 do i= 1, m
2309 denm=ts(i,j,k,is)/( 1.-rs(i,j,k,is)*rxa(i,j,k+1,1,is) )
2310 rra(i,j,k,1,is)=rr(i,j,k,is)+(td(i,j,k,is)
2311 * *rra(i,j,k+1,1,is)+tt(i,j,k,is)*rxa(i,j,k+1,1,is))*denm
2312 rxa(i,j,k,1,is)= rs(i,j,k,is)+ts(i,j,k,is)
2313 * *rxa(i,j,k+1,1,is)*denm
2314 rra(i,j,k,2,is)=rra(i,j,k,1,is)
2315 rxa(i,j,k,2,is)=rxa(i,j,k,1,is)
2316 enddo
2317 enddo
2318 enddo
2319
2320 c-----for middle clouds
2321
2322 do 20 im=1,2
2323
2324 do k= icb-1,ict,-1
2325 do j= 1, n
2326 do i= 1, m
2327 denm=ts(i,j,k,im)/( 1.-rs(i,j,k,im)*rxa(i,j,k+1,im,is) )
2328 rra(i,j,k,im,is)= rr(i,j,k,im)+(td(i,j,k,im)
2329 * *rra(i,j,k+1,im,is)+tt(i,j,k,im)*rxa(i,j,k+1,im,is))*denm
2330 rxa(i,j,k,im,is)= rs(i,j,k,im)+ts(i,j,k,im)
2331 * *rxa(i,j,k+1,im,is)*denm
2332 enddo
2333 enddo
2334 enddo
2335
2336 20 continue
2337
2338 c-----integration over eight sky situations.
2339 c ih, im, is denotes high, middle and low cloud groups.
2340
2341 do 100 ih=1,2
2342
2343 c-----clear portion
2344
2345 if(ih.eq.1) then
2346 do j=1,n
2347 do i=1,m
2348 ch(i,j)=1.0-cc(i,j,1)
2349 enddo
2350 enddo
2351
2352 else
2353
2354 c-----cloudy portion
2355
2356 do j=1,n
2357 do i=1,m
2358 ch(i,j)=cc(i,j,1)
2359 enddo
2360 enddo
2361
2362 endif
2363
2364 do 100 im=1,2
2365
2366 c-----clear portion
2367
2368 if(im.eq.1) then
2369
2370 do j=1,n
2371 do i=1,m
2372 cm(i,j)=ch(i,j)*(1.0-cc(i,j,2))
2373 enddo
2374 enddo
2375
2376 else
2377
2378 c-----cloudy portion
2379
2380 do j=1,n
2381 do i=1,m
2382 cm(i,j)=ch(i,j)*cc(i,j,2)
2383 enddo
2384 enddo
2385
2386 endif
2387
2388 do 100 is=1,2
2389
2390 c-----clear portion
2391
2392 if(is.eq.1) then
2393
2394 do j=1,n
2395 do i=1,m
2396 ct(i,j)=cm(i,j)*(1.0-cc(i,j,3))
2397 enddo
2398 enddo
2399
2400 else
2401
2402 c-----cloudy portion
2403
2404 do j=1,n
2405 do i=1,m
2406 ct(i,j)=cm(i,j)*cc(i,j,3)
2407 enddo
2408 enddo
2409
2410 endif
2411
2412 c-----add one layer at a time, going down.
2413
2414 do k= icb, np
2415 do j= 1, n
2416 do i= 1, m
2417 denm = ts(i,j,k,is)/( 1.-rsa(i,j,k-1,ih,im)*rs(i,j,k,is) )
2418 tda(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*td(i,j,k,is)
2419 tta(i,j,k,ih,im)= tda(i,j,k-1,ih,im)*tt(i,j,k,is)
2420 * +(tda(i,j,k-1,ih,im)*rr(i,j,k,is)
2421 * *rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2422 rsa(i,j,k,ih,im)= rs(i,j,k,is)+ts(i,j,k,is)
2423 * *rsa(i,j,k-1,ih,im)*denm
2424 enddo
2425 enddo
2426 enddo
2427
2428 c-----add one layer at a time, going up.
2429
2430 do k= ict-1,1,-1
2431 do j= 1, n
2432 do i= 1, m
2433 denm =ts(i,j,k,ih)/(1.-rs(i,j,k,ih)*rxa(i,j,k+1,im,is))
2434 rra(i,j,k,im,is)= rr(i,j,k,ih)+(td(i,j,k,ih)
2435 * *rra(i,j,k+1,im,is)+tt(i,j,k,ih)*rxa(i,j,k+1,im,is))*denm
2436 rxa(i,j,k,im,is)= rs(i,j,k,ih)+ts(i,j,k,ih)
2437 * *rxa(i,j,k+1,im,is)*denm
2438 enddo
2439 enddo
2440 enddo
2441
2442 c-----compute fluxes following eq (5) of Chou (1992)
2443
2444 c fdndir is the direct downward flux
2445 c fdndif is the diffuse downward flux
2446 c fupdif is the diffuse upward flux
2447
2448 do k=2,np+1
2449 do j=1, n
2450 do i=1, m
2451 denm= 1./(1.- rxa(i,j,k,im,is)*rsa(i,j,k-1,ih,im))
2452 fdndir(i,j)= tda(i,j,k-1,ih,im)
2453 xx = tda(i,j,k-1,ih,im)*rra(i,j,k,im,is)
2454 fdndif(i,j)= (xx*rsa(i,j,k-1,ih,im)+tta(i,j,k-1,ih,im))*denm
2455 fupdif= (xx+tta(i,j,k-1,ih,im)*rxa(i,j,k,im,is))*denm
2456 flxdn(i,j,k)=fdndir(i,j)+fdndif(i,j)-fupdif
2457 enddo
2458 enddo
2459 enddo
2460
2461 do j=1, n
2462 do i=1, m
2463 flxdn(i,j,1)=1.0-rra(i,j,1,im,is)
2464 enddo
2465 enddo
2466
2467 c-----summation of fluxes over all (eight) sky situations.
2468
2469 do k=1,np+1
2470 do j=1,n
2471 do i=1,m
2472 if(ih.eq.1 .and. im.eq.1 .and. is.eq.1) then
2473 fclr(i,j,k)=flxdn(i,j,k)
2474 endif
2475 fall(i,j,k)=fall(i,j,k)+flxdn(i,j,k)*ct(i,j)
2476 enddo
2477 enddo
2478 enddo
2479
2480 do j=1,n
2481 do i=1,m
2482 fsdir(i,j)=fsdir(i,j)+fdndir(i,j)*ct(i,j)
2483 fsdif(i,j)=fsdif(i,j)+fdndif(i,j)*ct(i,j)
2484 enddo
2485 enddo
2486
2487 100 continue
2488
2489 return
2490 end
2491
2492 c*****************************************************************
2493
2494 subroutine flxco2(m,n,np,swc,swh,csm,df)
2495
2496 c*****************************************************************
2497
2498 c-----compute the reduction of clear-sky downward solar flux
2499 c due to co2 absorption.
2500
2501 implicit none
2502
2503 c-----input parameters
2504
2505 integer m,n,np
2506 _RL csm(m,n),swc(m,n,np+1),swh(m,n,np+1),cah(22,19)
2507
2508 c-----output (undated) parameter
2509
2510 _RL df(m,n,np+1)
2511
2512 c-----temporary array
2513
2514 integer i,j,k,ic,iw
2515 _RL xx,clog1,wlog,dc,dw,x1,x2,y2
2516
2517 c********************************************************************
2518 c-----include co2 look-up table
2519
2520 #include "cah-dat.h"
2521 c save cah
2522
2523 c********************************************************************
2524 c-----table look-up for the reduction of clear-sky solar
2525 c radiation due to co2. The fraction 0.0343 is the
2526 c extraterrestrial solar flux in the co2 bands.
2527
2528 do k= 2, np+1
2529 do j= 1, n
2530 do i= 1, m
2531 xx=1./.3
2532 clog1=log10(swc(i,j,k)*csm(i,j))
2533 wlog=log10(swh(i,j,k)*csm(i,j))
2534 ic=int( (clog1+3.15)*xx+1.)
2535 iw=int( (wlog+4.15)*xx+1.)
2536 if(ic.lt.2)ic=2
2537 if(iw.lt.2)iw=2
2538 if(ic.gt.22)ic=22
2539 if(iw.gt.19)iw=19
2540 dc=clog1-float(ic-2)*.3+3.
2541 dw=wlog-float(iw-2)*.3+4.
2542 x1=cah(1,iw-1)+(cah(1,iw)-cah(1,iw-1))*xx*dw
2543 x2=cah(ic-1,iw-1)+(cah(ic-1,iw)-cah(ic-1,iw-1))*xx*dw
2544 y2=x2+(cah(ic,iw-1)-cah(ic-1,iw-1))*xx*dc
2545 if (x1.lt.y2) x1=y2
2546 df(i,j,k)=df(i,j,k)+0.0343*(x1-y2)
2547 enddo
2548 enddo
2549 enddo
2550
2551 return
2552 end

  ViewVC Help
Powered by ViewVC 1.1.22