/[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.13 - (show annotations) (download)
Mon Jul 26 19:51:08 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint54d_post
Changes since 1.12: +3 -2 lines
Debugging fizhi still

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

  ViewVC Help
Powered by ViewVC 1.1.22