/[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.14 - (show annotations) (download)
Wed Jul 28 01:25:07 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.13: +36 -8 lines
debugging

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

  ViewVC Help
Powered by ViewVC 1.1.22