/[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.24 - (show annotations) (download)
Thu Jun 16 16:46:12 2005 UTC (19 years ago) by ce107
Branch: MAIN
Changes since 1.23: +23 -23 lines
Fix fizhi constants to _d 0 form as the IBM XL compiler complains on mixed
precision calls to intrinsics like max, min. Only problematic cases have
been altered - consider compiling with -qdpc=e to fix the precision of the
rest.

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

  ViewVC Help
Powered by ViewVC 1.1.22