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

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

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


Revision 1.25 - (show annotations) (download)
Thu May 5 18:40:42 2005 UTC (19 years ago) by molod
Branch: MAIN
Changes since 1.24: +3 -3 lines
More cleaning - use number of land tiles to fill tcanopy - not total number of tiles

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/update_earth_exports.F,v 1.24 2005/03/09 00:35:09 molod Exp $
2 C $Name: $
3
4 #include "FIZHI_OPTIONS.h"
5 subroutine update_earth_exports (myTime, myIter, myThid)
6 c----------------------------------------------------------------------
7 c Subroutine update_earth_exports - 'Wrapper' routine to update
8 c the fields related to the earth's surface that are needed
9 c by fizhi.
10 c
11 c Call: getlgr (Set the leaf area index and surface greenness,
12 c based on veg type and month)
13 c getalb (Set the 4 albedos based on veg type, snow and time)
14 c getemiss (Set the surface emissivity based on the veg type
15 c and the snow depth)
16 c-----------------------------------------------------------------------
17 implicit none
18 #include "SIZE.h"
19 #include "GRID.h"
20 #include "fizhi_land_SIZE.h"
21 #include "fizhi_SIZE.h"
22 #include "fizhi_coms.h"
23 #include "chronos.h"
24 #include "gridalt_mapping.h"
25 #include "fizhi_land_coms.h"
26 #include "fizhi_earth_coms.h"
27 #include "fizhi_ocean_coms.h"
28 #include "EEPARAMS.h"
29
30 integer myIter, myThid
31 _RL myTime
32
33 logical alarm
34 external alarm
35 _RL lats(sNx,sNy), lons(sNx,sNy), cosz(sNx,sNy)
36 _RL fraci(sNx,sNy), fracl(sNx,sNy)
37 _RL ficetile(nchp)
38 _RL radius
39 _RL tmpij(sNx,sNy)
40 _RL tmpchp(nchp)
41 integer i, j, n, bi, bj
42 integer im1, im2, jm1, jm2, idim1, idim2, jdim1, jdim2
43 integer sec, day, month
44 integer nmonf,ndayf,nsecf
45 nsecf(n) = n/10000*3600 + mod(n,10000)/100* 60 + mod(n,100)
46 nmonf(n) = mod(n,10000)/100
47 ndayf(n) = mod(n,100)
48
49 idim1 = 1-OLx
50 idim2 = sNx+OLx
51 jdim1 = 1-OLy
52 jdim2 = sNy+OLy
53 im1 = 1
54 im2 = sNx
55 jm1 = 1
56 jm2 = sNy
57 month = nmonf(nymd)
58 day = ndayf(nymd)
59 sec = nsecf(nhms)
60
61 do bj = myByLo(myThid), myByHi(myThid)
62 do bi = myBxLo(myThid), myBxHi(myThid)
63 do j = jm1,jm2
64 do i = im1,im2
65 lons(i,j) = xC(i,j,bi,bj)
66 lats(i,j) = yC(i,j,bi,bj)
67 enddo
68 enddo
69
70 call get_landfrac(im2,jm2,nSx,nSy,bi,bj,maxtyp,surftype,tilefrac,
71 . fracl)
72
73 do j = jm1,jm2
74 do i = im1,im2
75 if(sice(i,j,bi,bj).gt.0.) then
76 fraci(i,j) = 1.
77 else
78 fraci(i,j) = 0.
79 endif
80 enddo
81 enddo
82
83 C***********************************************************************
84 C* Get Leaf-Area-Index and Greenness Index *
85 C***********************************************************************
86
87 if( alarm('turb') .or. alarm('radsw') ) then
88 call getlgr (sec,month,day,chlt,ityp,nchpland(bi,bj),
89 . nchp,nSx,nSy,bi,bj,alai,agrn )
90 endif
91
92 C **********************************************************************
93 C Compute Surface Albedo
94 C **********************************************************************
95
96 if( alarm('radsw') ) then
97 call astro(nymd,nhms,lats,lons,im2*jm2,cosz,radius)
98 call getalb(sec,month,day,cosz,snodep,fraci,fracl,im2,jm2,nchp,
99 . nchptot(bi,bj),nchpland(bi,bj),nSx,nSy,bi,bj,igrd,ityp,
100 . chfr,chlt,alai,agrn,
101 . albvisdr,albvisdf,albnirdr,albnirdf )
102 endif
103
104 C **********************************************************************
105 C Compute Surface Emissivity
106 C **********************************************************************
107
108 if( alarm('radlw') ) then
109 call grd2msc(fraci,im2,jm2,igrd,ficetile,nchp,nchptot(bi,bj))
110 call getemiss(fracl,im2,jm2,nchp,nchptot(bi,bj),nSx,nSy,bi,bj,
111 . igrd,ityp,chfr,snodep,ficetile,emiss)
112 endif
113
114 C*********************************************************************
115 C Ground Temperature Over Ocean is from SST array,
116 C Over land is from tcanopy
117 C*********************************************************************
118
119 do j = jm1,jm2
120 do i = im1,im2
121 tmpij(i,j) = 0.
122 enddo
123 enddo
124 do i = 1,nchpland(bi,bj)
125 tmpchp(i) = tcanopy(i,bi,bj)
126 enddo
127 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),tmpchp,
128 . nchp,nchpland(bi,bj),fracl,tmpij,im2,jm2)
129 do j = jm1,jm2
130 do i = im1,im2
131 tgz(i,j,bi,bj) = tmpij(i,j)
132 if(fracl(i,j).lt.0.3.and.sice(i,j,bi,bj).eq.0.0)
133 . tgz(i,j,bi,bj) = sst(i,j,bi,bj)
134 enddo
135 enddo
136
137 enddo
138 enddo
139
140 return
141 end
142
143 SUBROUTINE SIBALB ( AVISDR, ANIRDR, AVISDF, ANIRDF,
144 . VLAI, VGRN, ZTH, SNW, ITYP, IRUN )
145
146 C*********************************************************************
147 C The input list is as follows:
148 C VLAI: the leaf area index.
149 C VGRN: the greenness index.
150 C ZTH: The cosine of the solar zenith angle.
151 C SNW: Snow cover in meters water equivalent.
152 C ITYP: The surface type (grass, bare soil, etc.)
153 C IRUN: Number of tiles (same as used for SUBROUTINE TILE).
154 C
155 C The output list is as follows:
156 C
157 C AVISDR: visible, direct albedo.
158 C ANIRDR: near infra-red, direct albedo.
159 C AVISDF: visible, diffuse albedo.
160 C ANIRDF: near infra-red, diffuse albedo.
161 C*******************************************************************
162
163 IMPLICIT NONE
164
165 INTEGER IRUN
166 _RL AVISDR (IRUN), ANIRDR (IRUN), AVISDF (IRUN), ANIRDF (IRUN)
167 _RL VLAI(IRUN),VGRN (IRUN), SNW(IRUN)
168 _RL ZTH(IRUN)
169 INTEGER ITYP (IRUN)
170
171 _RL ALVDRS, ALIDRS
172 _RL ALVDRDL, ALIDRDL
173 _RL ALVDRDD, ALIDRDD
174 _RL ALVDRI, ALIDRI
175 _RL minval
176 external minval
177
178 C Albedo of soil for visible direct solar radiation.
179 PARAMETER ( ALVDRS = 0.100 )
180 C Albedo of soil for infra-red direct solar radiation.
181 PARAMETER ( ALIDRS = 0.200 )
182 C Albedo of light desert for visible direct solar radiation.
183 PARAMETER ( ALVDRDL = 0.300 )
184 C Albedo of light desert for infra-red direct solar radiation.
185 PARAMETER ( ALIDRDL = 0.350 )
186 C Albedo of dark desert for visible direct solar radiation.
187 PARAMETER ( ALVDRDD = 0.250 )
188 C Albedo of dark desert for infra-red direct solar radiation.
189 PARAMETER ( ALIDRDD = 0.300 )
190 C Albedo of ice for visible direct solar radiation.
191 PARAMETER ( ALVDRI = 0.800 )
192 C Albedo of ice for infra-red direct solar radiation.
193 PARAMETER ( ALIDRI = 0.800 )
194
195 * --------------------------------------------------------------------------------------------
196
197 INTEGER NTYPS
198 INTEGER NLAI
199 _RL ZERO, ONE
200 _RL EPSLN, BLAI, DLAI
201 _RL ALATRM
202 PARAMETER (NLAI = 14 )
203 PARAMETER (EPSLN = 1.E-6)
204 PARAMETER (BLAI = 0.5)
205 PARAMETER (DLAI = 0.5)
206 PARAMETER (ZERO=0., ONE=1.0)
207 PARAMETER (ALATRM = BLAI + (NLAI - 1) * DLAI - EPSLN)
208 PARAMETER (NTYPS=10)
209
210
211 C ITYP: Vegetation type as follows:
212 C 1: BROADLEAF EVERGREEN TREES
213 C 2: BROADLEAF DECIDUOUS TREES
214 C 3: NEEDLELEAF TREES
215 C 4: GROUND COVER
216 C 5: BROADLEAF SHRUBS
217 C 6: DWARF TREES (TUNDRA)
218 C 7: BARE SOIL
219 C 8: LIGHT DESERT
220 C 9: GLACIER
221 C 10: DARK DESERT
222 C
223
224 INTEGER I, LAI
225 _RL FAC,GAMMA,BETA,ALPHA,DX,DY,ALA,GRN (2),SNWALB(4,NTYPS)
226 _RL COEFF
227
228 _RL ALVDR (NLAI, 2, NTYPS)
229 _RL BTVDR (NLAI, 2, NTYPS)
230 _RL GMVDR (NLAI, 2, NTYPS)
231 _RL ALIDR (NLAI, 2, NTYPS)
232 _RL BTIDR (NLAI, 2, NTYPS)
233 _RL GMIDR (NLAI, 2, NTYPS)
234
235 C (Data statements for ALVDR described in full; data statements for
236 C other constants follow same framework.)
237
238 C BROADLEAF EVERGREEN (ITYP=4); GREEN=0.33; LAI: .5-7
239 DATA (ALVDR (I, 1, 1), I = 1, 14)
240 ` /0.0808, 0.0796, 0.0792, 0.0790, 10*0.0789/
241
242 C BROADLEAF EVERGREEN (ITYP=4); GREEN=0.67; LAI: .5-7
243 DATA (ALVDR (I, 2, 1), I = 1, 14)
244 ` /0.0788, 0.0775, 0.0771, 0.0769, 10*0.0768/
245
246 C BROADLEAF DECIDUOUS (ITYP=1); GREEN=0.33; LAI: .5-7
247 DATA (ALVDR (I, 1, 2), I = 1, 14)
248 ` /0.0803, 0.0790, 0.0785, 0.0784, 3*0.0783, 7*0.0782/
249
250 C BROADLEAF DECIDUOUS (ITYP=1); GREEN=0.67; LAI: .5-7
251 DATA (ALVDR (I, 2, 2), I = 1, 14)
252 ` /0.0782, 0.0770, 0.0765, 0.0763, 10*0.0762/
253
254 C NEEDLELEAF (ITYP=3); GREEN=0.33; LAI=.5-7
255 DATA (ALVDR (I, 1, 3), I = 1, 14)
256 ` /0.0758, 0.0746, 0.0742, 0.0740, 10*0.0739/
257
258 C NEEDLELEAF (ITYP=3); GREEN=0.67; LAI=.5-7
259 DATA (ALVDR (I, 2, 3), I = 1, 14)
260 ` /0.0683, 0.0672, 0.0667, 2*0.0665, 9*0.0664/
261
262 C GROUNDCOVER (ITYP=2); GREEN=0.33; LAI=.5-7
263 DATA (ALVDR (I, 1, 4), I = 1, 14)
264 ` /0.2436, 0.2470, 0.2486, 0.2494, 0.2498, 0.2500, 2*0.2501,
265 ` 6*0.2502
266 ` /
267 C GROUNDCOVER (ITYP=2); GREEN=0.67; LAI=.5-7
268 DATA (ALVDR (I, 2, 4), I = 1, 14) /14*0.1637/
269
270 C BROADLEAF SHRUBS (ITYP=5); GREEN=0.33,LAI=.5-7
271 DATA (ALVDR (I, 1, 5), I = 1, 14)
272 & /0.0807, 0.0798, 0.0794, 0.0792, 0.0792, 9*0.0791/
273
274 C BROADLEAF SHRUBS (ITYP=5); GREEN=0.67,LAI=.5-7
275 DATA (ALVDR (I, 2, 5), I = 1, 14)
276 & /0.0787, 0.0777, 0.0772, 0.0771, 10*0.0770/
277
278 C DWARF TREES, OR TUNDRA (ITYP=6); GREEN=0.33,LAI=.5-7
279 DATA (ALVDR (I, 1, 6), I = 1, 14)
280 & /0.0802, 0.0791, 0.0787, 0.0786, 10*0.0785/
281
282 C DWARF TREES, OR TUNDRA (ITYP=6); GREEN=0.67,LAI=.5-7
283 DATA (ALVDR (I, 2, 6), I = 1, 14)
284 & /0.0781, 0.0771, 0.0767, 0.0765, 0.0765, 9*0.0764/
285
286
287 C BARE SOIL
288 DATA (ALVDR (I, 1, 7), I = 1, 14) /14*ALVDRS/
289 DATA (ALVDR (I, 2, 7), I = 1, 14) /14*ALVDRS/
290
291 C LIGHT DESERT (SAHARA, EG)
292 DATA (ALVDR (I, 1, 8), I = 1, 14) /14*ALVDRDL/
293 DATA (ALVDR (I, 2, 8), I = 1, 14) /14*ALVDRDL/
294
295 C ICE
296 DATA (ALVDR (I, 1, 9), I = 1, 14) /14*ALVDRI/
297 DATA (ALVDR (I, 2, 9), I = 1, 14) /14*ALVDRI/
298
299 C DARK DESERT (AUSTRALIA, EG)
300 DATA (ALVDR (I, 1, 10), I = 1, 14) /14*ALVDRDD/
301 DATA (ALVDR (I, 2, 10), I = 1, 14) /14*ALVDRDD/
302 C****
303 C**** -------------------------------------------------
304 DATA (BTVDR (I, 1, 1), I = 1, 14)
305 ` /0.0153, 0.0372, 0.0506, 0.0587, 0.0630, 0.0652, 0.0663,
306 ` 0.0668, 0.0671, 0.0672, 4*0.0673
307 ` /
308 DATA (BTVDR (I, 2, 1), I = 1, 14)
309 * /0.0135, 0.0354, 0.0487, 0.0568, 0.0611, 0.0633, 0.0644,
310 ` 0.0650, 0.0652, 0.0654, 0.0654, 3*0.0655
311 ` /
312 DATA (BTVDR (I, 1, 2), I = 1, 14)
313 * /0.0148, 0.0357, 0.0462, 0.0524, 0.0554, 0.0569, 0.0576,
314 ` 0.0579, 0.0580, 0.0581, 0.0581, 3*0.0582
315 ` /
316 DATA (BTVDR (I, 2, 2), I = 1, 14)
317 * /0.0131, 0.0342, 0.0446, 0.0508, 0.0539, 0.0554, 0.0560,
318 ` 0.0564, 0.0565, 5*0.0566
319 ` /
320 DATA (BTVDR (I, 1, 3), I = 1, 14)
321 * /0.0108, 0.0334, 0.0478, 0.0571, 0.0624, 0.0652, 0.0666,
322 ` 0.0673, 0.0677, 0.0679, 4*0.0680
323 ` /
324 DATA (BTVDR (I, 2, 3), I = 1, 14)
325 * /0.0034, 0.0272, 0.0408, 0.0501, 0.0554, 0.0582, 0.0597,
326 * 0.0604, 0.0608, 0.0610, 4*0.0611
327 ` /
328 DATA (BTVDR (I, 1, 4), I = 1, 14)
329 * /0.2050, 0.2524, 0.2799, 0.2947, 0.3022, 0.3059, 0.3076,
330 * 0.3085, 0.3088, 0.3090, 4*0.3091
331 ` /
332 DATA (BTVDR (I, 2, 4), I = 1, 14)
333 * /0.1084, 0.1404, 0.1617, 0.1754, 0.1837, 0.1887, 0.1915,
334 * 0.1931, 0.1940, 0.1946, 0.1948, 0.1950, 2*0.1951
335 ` /
336 DATA (BTVDR (I, 1, 5), I = 1, 14)
337 & /0.0203, 0.0406, 0.0548, 0.0632, 0.0679, 0.0703, 0.0716,
338 & 0.0722, 0.0726, 0.0727, 0.0728, 0.0728, 0.0728, 0.0729
339 ` /
340
341 DATA (BTVDR (I, 2, 5), I = 1, 14)
342 & /0.0184, 0.0385, 0.0526, 0.0611, 0.0658, 0.0683, 0.0696,
343 & 0.0702, 0.0705, 0.0707, 4*0.0708
344 ` /
345
346 DATA (BTVDR (I, 1, 6), I = 1, 14)
347 & /0.0199, 0.0388, 0.0494, 0.0554, 0.0584, 0.0599, 0.0606,
348 & 0.0609, 0.0611, 5*0.0612
349 ` /
350
351 DATA (BTVDR (I, 2, 6), I = 1, 14)
352 & /0.0181, 0.0371, 0.0476, 0.0537, 0.0568, 0.0583, 0.0590,
353 & 0.0593, 0.0595, 0.0595, 4*0.0596
354 ` /
355
356 DATA (BTVDR (I, 1, 7), I = 1, 14) /14*0./
357 DATA (BTVDR (I, 2, 7), I = 1, 14) /14*0./
358
359 DATA (BTVDR (I, 1, 8), I = 1, 14) /14*0./
360 DATA (BTVDR (I, 2, 8), I = 1, 14) /14*0./
361
362 DATA (BTVDR (I, 1, 9), I = 1, 14) /14*0./
363 DATA (BTVDR (I, 2, 9), I = 1, 14) /14*0./
364
365 DATA (BTVDR (I, 1, 10), I = 1, 14) /14*0./
366 DATA (BTVDR (I, 2, 10), I = 1, 14) /14*0./
367
368 C****
369 C**** -----------------------------------------------------------
370 DATA (GMVDR (I, 1, 1), I = 1, 14)
371 ` /0.0814, 0.1361, 0.2078, 0.2650, 0.2986, 0.3169, 0.3265,
372 * 0.3313, 0.3337, 0.3348, 0.3354, 0.3357, 2*0.3358
373 ` /
374 DATA (GMVDR (I, 2, 1), I = 1, 14)
375 * /0.0760, 0.1336, 0.2034, 0.2622, 0.2969, 0.3159, 0.3259,
376 * 0.3309, 0.3333, 0.3346, 0.3352, 0.3354, 2*0.3356
377 ` /
378 DATA (GMVDR (I, 1, 2), I = 1, 14)
379 * /0.0834, 0.1252, 0.1558, 0.1927, 0.2131, 0.2237, 0.2290,
380 * 0.2315, 0.2327, 0.2332, 0.2335, 2*0.2336, 0.2337
381 ` /
382 DATA (GMVDR (I, 2, 2), I = 1, 14)
383 * /0.0789, 0.1235, 0.1531, 0.1912, 0.2122, 0.2232, 0.2286,
384 * 0.2312, 0.2324, 0.2330, 0.2333, 0.2334, 2*0.2335
385 ` /
386 DATA (GMVDR (I, 1, 3), I = 1, 14)
387 * /0.0647, 0.1342, 0.2215, 0.2968, 0.3432, 0.3696, 0.3838,
388 * 0.3912, 0.3950, 0.3968, 0.3978, 0.3982, 0.3984, 0.3985
389 ` /
390 DATA (GMVDR (I, 2, 3), I = 1, 14)
391 * /0.0258, 0.1227, 0.1999, 0.2825, 0.3339, 0.3634, 0.3794,
392 * 0.3877, 0.3919, 0.3940, 0.3950, 0.3956, 0.3958, 0.3959
393 ` /
394 DATA (GMVDR (I, 1, 4), I = 1, 14)
395 * /0.3371, 0.5762, 0.7159, 0.7927, 0.8324, 0.8526, 0.8624,
396 * 0.8671, 0.8693, 0.8704, 0.8709, 0.8710, 2*0.8712
397 ` /
398 DATA (GMVDR (I, 2, 4), I = 1, 14)
399 * /0.2634, 0.4375, 0.5532, 0.6291, 0.6763, 0.7048, 0.7213,
400 * 0.7310, 0.7363, 0.7395, 0.7411, 0.7420, 0.7426, 0.7428
401 ` /
402 DATA (GMVDR (I, 1, 5), I = 1, 14)
403 & /0.0971, 0.1544, 0.2511, 0.3157, 0.3548, 0.3768, 0.3886,
404 & 0.3948, 0.3978, 0.3994, 0.4001, 0.4006, 0.4007, 0.4008
405 ` /
406
407 DATA (GMVDR (I, 2, 5), I = 1, 14)
408 & /0.0924, 0.1470, 0.2458, 0.3123, 0.3527, 0.3756, 0.3877,
409 & 0.3942, 0.3974, 0.3990, 0.3998, 0.4002, 0.4004, 0.4005
410 ` /
411
412 DATA (GMVDR (I, 1, 6), I = 1, 14)
413 & /0.0970, 0.1355, 0.1841, 0.2230, 0.2447, 0.2561, 0.2617,
414 & 0.2645, 0.2658, 0.2664, 0.2667, 3*0.2669
415 ` /
416
417 DATA (GMVDR (I, 2, 6), I = 1, 14)
418 & /0.0934, 0.1337, 0.1812, 0.2213, 0.2437, 0.2554, 0.2613,
419 & 0.2642, 0.2656, 0.2662, 0.2665, 0.2667, 0.2667, 0.2668
420 ` /
421
422 DATA (GMVDR (I, 1, 7), I = 1, 14) /14*1./
423 DATA (GMVDR (I, 2, 7), I = 1, 14) /14*1./
424
425 DATA (GMVDR (I, 1, 8), I = 1, 14) /14*1./
426 DATA (GMVDR (I, 2, 8), I = 1, 14) /14*1./
427
428 DATA (GMVDR (I, 1, 9), I = 1, 14) /14*1./
429 DATA (GMVDR (I, 2, 9), I = 1, 14) /14*1./
430
431 DATA (GMVDR (I, 1, 10), I = 1, 14) /14*1./
432 DATA (GMVDR (I, 2, 10), I = 1, 14) /14*1./
433
434 C****
435 C**** -----------------------------------------------------------
436
437 DATA (ALIDR (I, 1, 1), I = 1, 14)
438 * /0.2867, 0.2840, 0.2828, 0.2822, 0.2819, 0.2818, 2*0.2817,
439 * 6*0.2816
440 ` /
441 DATA (ALIDR (I, 2, 1), I = 1, 14)
442 * /0.3564, 0.3573, 0.3577, 0.3580, 2*0.3581, 8*0.3582/
443 DATA (ALIDR (I, 1, 2), I = 1, 14)
444 * /0.2848, 0.2819, 0.2804, 0.2798, 0.2795, 2*0.2793, 7*0.2792/
445 DATA (ALIDR (I, 2, 2), I = 1, 14)
446 * /0.3544, 0.3550, 0.3553, 2*0.3555, 9*0.3556/
447 DATA (ALIDR (I, 1, 3), I = 1, 14)
448 * /0.2350, 0.2311, 0.2293, 0.2285, 0.2281, 0.2280, 8*0.2279/
449 DATA (ALIDR (I, 2, 3), I = 1, 14)
450 * /0.2474, 0.2436, 0.2418, 0.2410, 0.2406, 0.2405, 3*0.2404,
451 * 5*0.2403
452 ` /
453 DATA (ALIDR (I, 1, 4), I = 1, 14)
454 * /0.5816, 0.6157, 0.6391, 0.6556, 0.6673, 0.6758, 0.6820,
455 * 0.6866, 0.6899, 0.6924, 0.6943, 0.6956, 0.6966, 0.6974
456 ` /
457 DATA (ALIDR (I, 2, 4), I = 1, 14)
458 * /0.5489, 0.5770, 0.5955, 0.6079, 0.6163, 0.6221, 0.6261,
459 * 0.6288, 0.6308, 0.6321, 0.6330, 0.6337, 0.6341, 0.6344
460 ` /
461 DATA (ALIDR (I, 1, 5), I = 1, 14)
462 & /0.2845, 0.2837, 0.2832, 0.2831, 0.2830, 9*0.2829/
463 DATA (ALIDR (I, 2, 5), I = 1, 14)
464 & /0.3532, 0.3562, 0.3578, 0.3586, 0.3590, 0.3592, 0.3594,
465 & 0.3594, 0.3594, 5*0.3595
466 ` /
467 DATA (ALIDR (I, 1, 6), I = 1, 14)
468 & /0.2825, 0.2812, 0.2806, 0.2803, 0.2802, 9*0.2801/
469 DATA (ALIDR (I, 2, 6), I = 1, 14)
470 & /0.3512, 0.3538, 0.3552, 0.3559, 0.3562, 0.3564, 0.3565,
471 & 0.3565, 6*0.3566
472 ` /
473
474 DATA (ALIDR (I, 1, 7), I = 1, 14) /14*ALIDRS/
475 DATA (ALIDR (I, 2, 7), I = 1, 14) /14*ALIDRS/
476
477 DATA (ALIDR (I, 1, 8), I = 1, 14) /14*ALIDRDL/
478 DATA (ALIDR (I, 2, 8), I = 1, 14) /14*ALIDRDL/
479
480 DATA (ALIDR (I, 1, 9), I = 1, 14) /14*ALIDRI/
481 DATA (ALIDR (I, 2, 9), I = 1, 14) /14*ALIDRI/
482
483 DATA (ALIDR (I, 1, 10), I = 1, 14) /14*ALIDRDD/
484 DATA (ALIDR (I, 2, 10), I = 1, 14) /14*ALIDRDD/
485
486 C****
487 C**** -----------------------------------------------------------
488 DATA (BTIDR (I, 1, 1), I = 1, 14)
489 * /0.1291, 0.1707, 0.1969, 0.2125, 0.2216, 0.2267, 0.2295,
490 * 0.2311, 0.2319, 0.2323, 0.2326, 2*0.2327, 0.2328
491 ` /
492 DATA (BTIDR (I, 2, 1), I = 1, 14)
493 * /0.1939, 0.2357, 0.2598, 0.2735, 0.2810, 0.2851, 0.2874,
494 * 0.2885, 0.2892, 0.2895, 0.2897, 3*0.2898
495 ` /
496 DATA (BTIDR (I, 1, 2), I = 1, 14)
497 * /0.1217, 0.1522, 0.1713, 0.1820, 0.1879, 0.1910, 0.1926,
498 * 0.1935, 0.1939, 0.1942, 2*0.1943, 2*0.1944
499 ` /
500 DATA (BTIDR (I, 2, 2), I = 1, 14)
501 * /0.1781, 0.2067, 0.2221, 0.2301, 0.2342, 0.2363, 0.2374,
502 * 0.2379, 0.2382, 0.2383, 2*0.2384, 2*0.2385
503 ` /
504 DATA (BTIDR (I, 1, 3), I = 1, 14)
505 * /0.0846, 0.1299, 0.1614, 0.1814, 0.1935, 0.2004, 0.2043,
506 * 0.2064, 0.2076, 0.2082, 0.2085, 2*0.2087, 0.2088
507 ` /
508 DATA (BTIDR (I, 2, 3), I = 1, 14)
509 * /0.0950, 0.1410, 0.1722, 0.1921, 0.2042, 0.2111, 0.2151,
510 * 0.2172, 0.2184, 0.2191, 0.2194, 0.2196, 2*0.2197
511 ` /
512 DATA (BTIDR (I, 1, 4), I = 1, 14)
513 * /0.5256, 0.7444, 0.9908, 1.2700, 1.5680, 1.8505, 2.0767,
514 * 2.2211, 2.2808, 2.2774, 2.2362, 2.1779, 2.1160, 2.0564
515 ` /
516 DATA (BTIDR (I, 2, 4), I = 1, 14)
517 * /0.4843, 0.6714, 0.8577, 1.0335, 1.1812, 1.2858, 1.3458,
518 * 1.3688, 1.3685, 1.3546, 1.3360, 1.3168, 1.2989, 1.2838
519 ` /
520 DATA (BTIDR (I, 1, 5), I = 1, 14)
521 & /0.1498, 0.1930, 0.2201, 0.2364, 0.2460, 0.2514, 0.2544,
522 & 0.2560, 0.2569, 0.2574, 0.2577, 0.2578, 0.2579, 0.2579
523 ` /
524
525 DATA (BTIDR (I, 2, 5), I = 1, 14)
526 & /0.2184, 0.2656, 0.2927, 0.3078, 0.3159, 0.3202, 0.3224,
527 & 0.3235, 0.3241, 0.3244, 0.3245, 3*0.3246
528 ` /
529
530 DATA (BTIDR (I, 1, 6), I = 1, 14)
531 & /0.1369, 0.1681, 0.1860, 0.1958, 0.2010, 0.2038, 0.2053,
532 & 0.2060, 0.2064, 0.2066, 0.2067, 3*0.2068
533 ` /
534
535 DATA (BTIDR (I, 2, 6), I = 1, 14)
536 & /0.1969, 0.2268, 0.2416, 0.2488, 0.2521, 0.2537, 0.2544,
537 & 0.2547, 0.2548, 5*0.2549
538 ` /
539
540
541 DATA (BTIDR (I, 1, 7), I = 1, 14) /14*0./
542 DATA (BTIDR (I, 2, 7), I = 1, 14) /14*0./
543
544 DATA (BTIDR (I, 1, 8), I = 1, 14) /14*0./
545 DATA (BTIDR (I, 2, 8), I = 1, 14) /14*0./
546
547 DATA (BTIDR (I, 1, 9), I = 1, 14) /14*0./
548 DATA (BTIDR (I, 2, 9), I = 1, 14) /14*0./
549
550 DATA (BTIDR (I, 1, 10), I = 1, 14) /14*0./
551 DATA (BTIDR (I, 2, 10), I = 1, 14) /14*0./
552
553 C****
554 C**** --------------------------------------------------------------
555 DATA (GMIDR (I, 1, 1), I = 1, 14)
556 * /0.1582, 0.2581, 0.3227, 0.3635, 0.3882, 0.4026, 0.4108,
557 * 0.4154, 0.4179, 0.4193, 0.4200, 0.4204, 0.4206, 0.4207
558 ` /
559 DATA (GMIDR (I, 2, 1), I = 1, 14)
560 * /0.1934, 0.3141, 0.3818, 0.4200, 0.4415, 0.4533, 0.4598,
561 * 0.4633, 0.4651, 0.4662, 0.4667, 0.4671, 2*0.4672
562 ` /
563 DATA (GMIDR (I, 1, 2), I = 1, 14)
564 * /0.1347, 0.1871, 0.2277, 0.2515, 0.2651, 0.2727, 0.2768,
565 * 0.2790, 0.2801, 0.2808, 0.2811, 0.2812, 0.2813, 0.2814
566 ` /
567 DATA (GMIDR (I, 2, 2), I = 1, 14)
568 * /0.1440, 0.2217, 0.2629, 0.2839, 0.2947, 0.3003, 0.3031,
569 * 0.3046, 0.3054, 0.3058, 0.3060, 2*0.3061, 0.3062
570 ` /
571 DATA (GMIDR (I, 1, 3), I = 1, 14)
572 * /0.1372, 0.2368, 0.3235, 0.3839, 0.4229, 0.4465, 0.4602,
573 * 0.4679, 0.4722, 0.4745, 0.4758, 0.4764, 0.4768, 0.4770
574 ` /
575 DATA (GMIDR (I, 2, 3), I = 1, 14)
576 * /0.1435, 0.2524, 0.3370, 0.3955, 0.4332, 0.4563, 0.4697,
577 * 0.4773, 0.4815, 0.4839, 0.4851, 0.4858, 0.4861, 0.4863
578 ` /
579 DATA (GMIDR (I, 1, 4), I = 1, 14)
580 * /0.4298, 0.9651, 1.6189, 2.4084, 3.2992, 4.1928, 4.9611,
581 * 5.5095, 5.8085, 5.9069, 5.8726, 5.7674, 5.6346, 5.4944
582 ` /
583 DATA (GMIDR (I, 2, 4), I = 1, 14)
584 * /0.4167, 0.8974, 1.4160, 1.9414, 2.4147, 2.7803, 3.0202,
585 * 3.1468, 3.1954, 3.1932, 3.1676, 3.1328, 3.0958, 3.0625
586 ` /
587 DATA (GMIDR (I, 1, 5), I = 1, 14)
588 & /0.1959, 0.3203, 0.3985, 0.4472, 0.4766, 0.4937, 0.5034,
589 & 0.5088, 0.5117, 0.5134, 0.5143, 0.5147, 0.5150, 0.5152
590 ` /
591
592 DATA (GMIDR (I, 2, 5), I = 1, 14)
593 & /0.2328, 0.3859, 0.4734, 0.5227, 0.5498, 0.5644, 0.5720,
594 & 0.5761, 0.5781, 0.5792, 0.5797, 0.5800, 0.5802, 0.5802
595 ` /
596
597 DATA (GMIDR (I, 1, 6), I = 1, 14)
598 & /0.1447, 0.2244, 0.2698, 0.2953, 0.3094, 0.3170, 0.3211,
599 & 0.3233, 0.3244, 0.3250, 0.3253, 0.3255, 0.3256, 0.3256
600 ` /
601
602 DATA (GMIDR (I, 2, 6), I = 1, 14)
603 & /0.1643, 0.2624, 0.3110, 0.3347, 0.3461, 0.3517, 0.3543,
604 & 0.3556, 0.3562, 0.3564, 0.3565, 0.3566, 0.3566, 0.3566
605 ` /
606
607 DATA (GMIDR (I, 1, 7), I = 1, 14) /14*1./
608 DATA (GMIDR (I, 2, 7), I = 1, 14) /14*1./
609
610 DATA (GMIDR (I, 1, 8), I = 1, 14) /14*1./
611 DATA (GMIDR (I, 2, 8), I = 1, 14) /14*1./
612
613 DATA (GMIDR (I, 1, 9), I = 1, 14) /14*1./
614 DATA (GMIDR (I, 2, 9), I = 1, 14) /14*1./
615
616 DATA (GMIDR (I, 1, 10), I = 1, 14) /14*1./
617 DATA (GMIDR (I, 2, 10), I = 1, 14) /14*1./
618
619
620 C**** -----------------------------------------------------------
621
622 DATA GRN /0.33, 0.67/
623
624 #include "snwmid.h"
625 DATA SNWALB /.65, .38, .65, .38,
626 * .65, .38, .65, .38,
627 * .65, .38, .65, .38,
628 * .65, .38, .65, .38,
629 * .65, .38, .65, .38,
630 & .65, .38, .65, .38,
631 & .65, .38, .65, .38,
632 & .65, .38, .65, .38,
633 & .80, .60, .80, .60,
634 & .65, .38, .65, .38
635 ` /
636
637 #ifdef CRAY
638 #ifdef f77
639 cfpp$ expand (coeff)
640 #endif
641 #endif
642
643 DO 100 I=1,IRUN
644 ALA = MIN (MAX (ZERO, VLAI(I)), ALATRM)
645 LAI = 1 + MAX(0, INT((ALA-BLAI)/DLAI) )
646 DX = (ALA - (BLAI+(LAI-1)*DLAI)) * (ONE/DLAI)
647 DY = (VGRN(I)- GRN(1)) * (ONE/(GRN(2) - GRN(1)))
648
649 ALPHA = COEFF (ALVDR (1, 1, ITYP (I)), NLAI, LAI ,DX, DY)
650 BETA = COEFF (BTVDR (1, 1, ITYP (I)), NLAI, LAI ,DX, DY)
651 GAMMA = COEFF (GMVDR (1, 1, ITYP (I)), NLAI, LAI ,DX, DY)
652
653 AVISDR(I) = ALPHA - ZTH(I)*BETA / (GAMMA+ZTH(I))
654 AVISDF(I) = ALPHA-BETA
655 * + 2.*BETA*GAMMA*(1.-GAMMA*LOG((1.+GAMMA)/GAMMA))
656
657 ALPHA = COEFF (ALIDR (1, 1, ITYP (I)), NLAI, LAI ,DX, DY)
658 BETA = COEFF (BTIDR (1, 1, ITYP (I)), NLAI, LAI ,DX, DY)
659 GAMMA = COEFF (GMIDR (1, 1, ITYP (I)), NLAI, LAI ,DX, DY)
660
661 ANIRDR(I) = ALPHA - ZTH(I)*BETA / (GAMMA+ZTH(I))
662 ANIRDF(I) = ALPHA-BETA
663 * + 2.*BETA*GAMMA*(1.-GAMMA*LOG((1.+GAMMA)/GAMMA))
664
665 IF (SNW (I) .GT. ZERO) THEN
666 FAC = SNW(I) / (SNW(I) + SNWMID(ITYP(I)))
667
668 AVISDR(I) = AVISDR(I) + (SNWALB(1,ITYP(I)) - AVISDR(I)) * FAC
669 ANIRDR(I) = ANIRDR(I) + (SNWALB(2,ITYP(I)) - ANIRDR(I)) * FAC
670 AVISDF(I) = AVISDF(I) + (SNWALB(3,ITYP(I)) - AVISDF(I)) * FAC
671 ANIRDF(I) = ANIRDF(I) + (SNWALB(4,ITYP(I)) - ANIRDF(I)) * FAC
672 ENDIF
673
674 100 CONTINUE
675
676 RETURN
677 END
678 FUNCTION COEFF(TABLE, NTABL, LAI ,DX, DY)
679
680 INTEGER NTABL, LAI
681 _RL coeff
682 _RL TABLE (NTABL, 2), DX, DY
683
684 COEFF = (TABLE(LAI, 1)
685 * + (TABLE(LAI ,2) - TABLE(LAI ,1)) * DY ) * (1.0-DX)
686 * + (TABLE(LAI+1,1)
687 * + (TABLE(LAI+1,2) - TABLE(LAI+1,1)) * DY ) * DX
688
689 RETURN
690 END
691
692 SUBROUTINE GETLGR(sec,IMON,IDAY,ALAT,ITYP,NCHPS,nchpdim,
693 . nSx,nSy,bi,bj,ALAI,AGRN)
694 C*********************************************************************
695 implicit none
696
697 integer ntyps
698 _RL one,daylen
699 PARAMETER (NTYPS=10)
700 parameter (one = 1.)
701 parameter (daylen = 86400.)
702
703 integer sec, imon, iday, nchps, nchpdim, nSx, nSy, bi, bj
704 _RL ALAI(nchpdim,nSx,nSy), AGRN(nchpdim,nSx,nSy)
705 _RL ALAT(nchpdim)
706 integer ITYP(nchpdim,nSx,nSy)
707
708 integer i,midmon,midm,midp,id,k1,k2,kk1,kk2
709 _RL fac
710
711 INTEGER DAYS(12)
712 DATA DAYS/31,28,31,30,31,30,31,31,30,31,30,31/
713
714 _RL VGLA(12,NTYPS), VGGR(12,NTYPS)
715
716 DATA VGLA /
717 1 5.117, 5.117, 5.117, 5.117, 5.117, 5.117, 5.117, 5.117,
718 1 5.117, 5.117, 5.117, 5.117,
719 2 0.520, 0.520, 0.867, 2.107, 4.507, 6.773, 7.173, 6.507,
720 2 5.040, 2.173, 0.867, 0.520,
721 3 8.760, 9.160, 9.827,10.093,10.360,10.760,10.493,10.227,
722 3 10.093, 9.827, 9.160, 8.760,
723 4 0.782, 0.893, 1.004, 1.116, 1.782, 3.671, 4.782, 4.227,
724 4 2.004, 1.227, 1.004, 0.893,
725 5 3.760, 3.760, 2.760, 1.760, 1.760, 1.760, 1.760, 5.760,
726 5 10.760, 7.760, 4.760, 3.760,
727 6 0.739, 0.739, 0.739, 0.739, 0.739, 1.072, 5.072, 5.739,
728 6 4.405, 0.739, 0.739, 0.739,
729 7 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
730 7 0.001, 0.001, 0.001, 0.001,
731 8 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
732 8 0.001, 0.001, 0.001, 0.001,
733 9 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
734 9 0.001, 0.001, 0.001, 0.001,
735 1 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
736 1 0.001, 0.001, 0.001, 0.001
737 & /
738
739
740 DATA VGGR
741 1 /0.905, 0.905, 0.905, 0.905, 0.905, 0.905, 0.905, 0.905,
742 1 0.905, 0.905, 0.905, 0.905,
743 2 0.026, 0.026, 0.415, 0.759, 0.888, 0.925, 0.836, 0.697,
744 2 0.331, 0.166, 0.015, 0.026,
745 3 0.913, 0.917, 0.923, 0.925, 0.927, 0.905, 0.902, 0.913,
746 3 0.898, 0.855, 0.873, 0.913,
747 4 0.568, 0.622, 0.664, 0.697, 0.810, 0.908, 0.813, 0.394,
748 4 0.443, 0.543, 0.553, 0.498,
749 5 0.798, 0.532, 0.362, 0.568, 0.568, 0.568, 0.568, 0.868,
750 5 0.651, 0.515, 0.630, 0.798,
751 6 0.451, 0.451, 0.451, 0.451, 0.451, 0.622, 0.920, 0.697,
752 6 0.076, 0.451, 0.451, 0.451,
753 7 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
754 7 0.001, 0.001, 0.001, 0.001,
755 8 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
756 8 0.001, 0.001, 0.001, 0.001,
757 9 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
758 9 0.001, 0.001, 0.001, 0.001,
759 1 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001, 0.001,
760 1 0.001, 0.001, 0.001, 0.001
761 & /
762
763
764 MIDMON = DAYS(IMON)/2 + 1
765
766
767 IF (IDAY .LT. MIDMON) THEN
768 K2 = IMON
769 K1 = MOD(IMON+10,12) + 1
770 ELSE
771 K1 = IMON
772 K2 = MOD(IMON,12) + 1
773 ENDIF
774
775 IF (IDAY .LT. MIDMON) THEN
776 MIDM = DAYS(K1)/2 + 1
777 MIDP = DAYS(K1) + MIDMON
778 ID = IDAY + DAYS(K1)
779 ELSE
780 MIDM = MIDMON
781 MIDP = DAYS(K2)/2 + 1 + DAYS(K1)
782 ID = IDAY
783 ENDIF
784
785 FAC = (float(ID -MIDM)*DAYLEN + SEC) /
786 * (float(MIDP-MIDM)*DAYLEN )
787
788 DO 220 I=1,NCHPS
789
790 IF(ALAT(I).GT.0.) THEN
791 KK1 = K1
792 KK2 = K2
793 ELSE
794 KK1 = MOD(K1+5,12) + 1
795 KK2 = MOD(K2+5,12) + 1
796 ENDIF
797
798 ALAI(I,bi,bj) = VGLA(KK2,ITYP(I,bi,bj))*FAC+
799 . VGLA(KK1,ITYP(I,bi,bj))*(ONE-FAC)
800 AGRN(I,bi,bj) = VGGR(KK2,ITYP(I,bi,bj))*FAC+
801 . VGGR(KK1,ITYP(I,bi,bj))*(ONE-FAC)
802
803 220 CONTINUE
804
805 RETURN
806 END
807
808 subroutine getalb(sec,month,day,cosz,snodep,fraci,fracg,im,jm,
809 . nchp,nchptot,nchpland,nSx,nSy,bi,bj,igrd,ityp,chfr,chlt,
810 . alai,agrn,albvr,albvf,albnr,albnf)
811 C***********************************************************************
812 C PURPOSE
813 C To act as an interface to routine sibalb, which calculates
814 C the four albedos for use by the shortwave radiation routine
815 C
816 C INPUT:
817 C sec - number of seconds into the day of current time
818 C month - month of the year of current time
819 C day - day of the month of current time
820 C cosz - local cosine of the zenith angle [im,jm]
821 C snodep - snow cover in meters [nchp,nSx,nSy]
822 C fraci - real array in grid space of total sea ice fraction [im,jm]
823 C fracg - real array in grid space of total land fraction [im,jm]
824 C im - model grid longitude dimension
825 C jm - model grid latitude dimension (number of lat. points)
826 C nchp - integer actual number of tiles in tile space
827 C nchpland - integer number of land tiles
828 C nSx - number of processors in x-direction
829 C nSy - number of processors in y-direction
830 C bi - processors index in x-direction
831 C bj - processors index in y-direction
832 C igrd - integer array in tile space of grid point number for each
833 C tile [nchp,nSx,nSy]
834 C ityp - integer array in tile space of land surface type for each
835 C tile [nchp,nSx,nSy]
836 C chfr - real array in tile space of land surface type fraction for
837 C each tile [nchp,nSx,nSy]
838 C chlt - real array in tile space of latitude value for each tile
839 C [nchp,nSx,nSy]
840 C
841 C OUTPUT:
842 C albvr - real array [im,jm] of visible direct beam albedo
843 C albvf - real array [im,jm] of visible diffuse beam albedo
844 C albnr - real array [im,jm] of near-ir direct beam albedo
845 C albnf - real array [im,jm] of near-ir diffuse beam albedo
846 C
847 C***********************************************************************
848 implicit none
849
850 integer sec,month,day,im,jm,nchp,nchptot,nchpland,nSx,nSy,bi,bj
851 _RL cosz(im,jm),fraci(im,jm),fracg(im,jm)
852 _RL snodep(nchp,nSx,nSy),chfr(nchp,nSx,nSy),chlt(nchp,nSx,nSy)
853 integer igrd(nchp,nSx,nSy),ityp(nchp,nSx,nSy)
854 _RL alai(nchp,nSx,nSy),agrn(nchp,nSx,nSy)
855 _RL albvr(im,jm,nSx,nSy),albvf(im,jm,nSx,nSy)
856 _RL albnr(im,jm,nSx,nSy),albnf(im,jm,nSx,nSy)
857
858 _RL one,a0,a1,a2,a3,ocnalb,albsi
859 PARAMETER (one = 1.)
860 PARAMETER (A0= 0.40670980)
861 PARAMETER (A1=-1.2523634 )
862 PARAMETER (A2= 1.4224051 )
863 PARAMETER (A3=-0.55573341)
864 PARAMETER (OCNALB=0.08)
865 PARAMETER (ALBSI=0.7)
866
867 _RL alboc(im,jm)
868 _RL AVISDR(nchp),ANIRDR(nchp),AVISDF(nchp)
869 _RL ANIRDF(nchp)
870 _RL zenith(nchp)
871 _RL tmpij(im,jm)
872 integer i,j
873
874 DO I=1,IM
875 DO J=1,JM
876 ALBOC(I,J) = A0 + (A1 + (A2 + A3*cosz(I,J))*cosz(I,J))*cosz(I,J)
877 ALBVR(I,J,bi,bj) = ALBSI*FRACI(I,J) + ALBOC(I,J)*(ONE-FRACI(I,J))
878 ALBNR(I,J,bi,bj) = ALBVR(I,J,bi,bj)
879 ALBVF(I,J,bi,bj) = ALBSI * FRACI(I,J) + OCNALB * (ONE-FRACI(I,J))
880 ALBNF(I,J,bi,bj) = ALBVF(I,J,bi,bj)
881 ENDDO
882 ENDDO
883
884 C and now some conversions from grid space to tile space before sibalb
885
886 call grd2msc(cosz,im,jm,igrd,zenith,nchp,nchpland)
887
888 C and now call sibalb
889
890 call sibalb(avisdr,anirdr,avisdf,anirdf,alai(1,bi,bj),
891 . agrn(1,bi,bj),zenith,snodep(1,bi,bj),ityp(1,bi,bj),nchpland)
892
893 C finally some transformations back to grid space for albedos
894
895 DO I=1,IM
896 DO J=1,JM
897 tmpij(i,j) = albvr(i,j,bi,bj)
898 ENDDO
899 ENDDO
900 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),avisdr,nchp,nchpland,
901 . fracg,tmpij,im,jm)
902
903 DO I=1,IM
904 DO J=1,JM
905 albvr(i,j,bi,bj) = tmpij(i,j)
906 ENDDO
907 ENDDO
908 DO I=1,IM
909 DO J=1,JM
910 tmpij(i,j) = albvf(i,j,bi,bj)
911 ENDDO
912 ENDDO
913 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),avisdf,nchp,nchpland,
914 . fracg,tmpij,im,jm)
915 DO I=1,IM
916 DO J=1,JM
917 albvf(i,j,bi,bj) = tmpij(i,j)
918 ENDDO
919 ENDDO
920 DO I=1,IM
921 DO J=1,JM
922 tmpij(i,j) = albnr(i,j,bi,bj)
923 ENDDO
924 ENDDO
925 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),anirdr,nchp,nchpland,
926 . fracg,tmpij,im,jm)
927 DO I=1,IM
928 DO J=1,JM
929 albnr(i,j,bi,bj) = tmpij(i,j)
930 ENDDO
931 ENDDO
932 DO I=1,IM
933 DO J=1,JM
934 tmpij(i,j) = albnf(i,j,bi,bj)
935 ENDDO
936 ENDDO
937 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),anirdf,nchp,nchpland,
938 . fracg,tmpij,im,jm)
939 DO I=1,IM
940 DO J=1,JM
941 albnf(i,j,bi,bj) = tmpij(i,j)
942 ENDDO
943 ENDDO
944
945 return
946 end
947
948 subroutine getemiss(fracg,im,jm,nchp,nchptot,nSx,nSy,bi,bj,
949 . igrd,ityp,chfr,snowdep,fraci,emiss)
950 C***********************************************************************
951 C PURPOSE
952 C To act as an interface to routine to emissivity, which calculates
953 C ten bands of surface emissivities for use by the longwave radiation
954 C
955 C INPUT:
956 C fracg - real array in grid space of total land fraction [im,jm]
957 C im - model grid longitude dimension
958 C jm - model grid latitude dimension (number of lat. points)
959 C nchp - integer actual number of tiles in tile space
960 C nSx - number of processors in x-direction
961 C nSy - number of processors in y-direction
962 C bi - processors index in x-direction
963 C bj - processors index in y-direction
964 C igrd - integer array in tile space of grid point number for each
965 C tile [nchp]
966 C ityp - integer array in tile space of land surface type for each
967 C tile [nchp]
968 C chfr - real array in tile space of land surface type fraction for
969 C each tile [nchp]
970 C snowdep - real array in tile space of snow depth (liquid water equiv)
971 C in mm [nchp]
972 C fraci - real array in tile space of sea ice fraction [nchp]
973 C
974 C OUTPUT:
975 C emiss - real array [im,jm,10,nSx,nSy] - surface emissivity (frac)
976 C
977 C***********************************************************************
978 implicit none
979 integer im,jm,nchp,nchptot,nSx,nSy,bi,bj
980 _RL fracg(im,jm)
981 _RL chfr(nchp,nSx,nSy)
982 integer igrd(nchp,nSx,nSy), ityp(nchp,nSx,nSy)
983 _RL snowdep(nchp,nSx,nSy)
984 _RL fraci(nchp)
985 _RL emiss(im,jm,10,nSx,nSy)
986
987 _RL emisstile(nchp,10)
988 _RL tmpij(im,jm)
989 integer i,j,k,n
990
991 do i = 1,10
992 do n = 1,nchptot
993 emisstile(n,i) = 1.
994 enddo
995 enddo
996
997 c call emissivity to get values in tile space
998 c -------------------------------------------
999 call emissivity(snowdep(1,bi,bj),fraci,nchp,nchptot,ityp(1,bi,bj),
1000 . emisstile)
1001
1002 c transform back to grid space for emissivities
1003 c ---------------------------------------------
1004 do k = 1,10
1005 do j = 1,jm
1006 do i = 1,im
1007 tmpij(i,j) = 0.0
1008 enddo
1009 enddo
1010 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),emisstile(1,k),nchp,
1011 . nchptot,fracg,tmpij,im,jm)
1012 do j = 1,jm
1013 do i = 1,im
1014 emiss(i,j,k,bi,bj) = tmpij(i,j)
1015 enddo
1016 enddo
1017 enddo
1018
1019 return
1020 end
1021
1022 subroutine emissivity (snowdepth,fraci,nchp,numpts,ityp,newemis)
1023 implicit none
1024 integer nchp,numpts
1025 integer ityp(nchp)
1026 _RL snowdepth(nchp)
1027 _RL fraci(nchp)
1028 _RL newemis(nchp,10)
1029
1030 _RL emis(12,11)
1031 _RL fac
1032 integer i,j
1033
1034 c-----------------------------------------------------------------------
1035 c NOTE: Emissivities were obtained for the following surface types:
1036 c ( 1) evergreen needleleaf = conifer
1037 c ( 2) evergreen broadleaf = conifer
1038 c ( 3) deciduous needleleaf = deciduous
1039 c ( 4) deciduous broadleaf = deciduous
1040 c ( 5) mixed forests = 1/2 conifer + 1/2 deciduous = tree
1041 c ( 6) closed shrublands = 3/4 tree + 1/4 quartz
1042 c ( 7) open shrubland = 1/4 tree + 3/4 quartz
1043 c ( 8) woody savannas = grass
1044 c ( 9) savannas = grass
1045 c (10) grasslands = grass
1046 c (11) permanent wetlands = 1/2 grass + 1/2 water
1047 c (12) croplands = grass
1048 c (13) urban = black body
1049 c (14) mosaic = 1/2 grass + 1/2 mixed forest
1050 c (15) snow/ice
1051 c (16) barren/sparsely vegetated = desert(quartz)
1052 c (17) water
1053 c (18) tundra = frost
1054 c
1055 c NOTE: Translation to Koster-Suarez surface types was as follows:
1056 c ( 1) broadleaf evergreen FROM above type 1 (conifer)
1057 c ( 2) broadleaf deciduous FROM above type 2 (deciduous)
1058 c ( 3) needleleaf evergreen FROM above type 1 (conifer)
1059 c ( 4) groundcover FROM above type 10 (grass)
1060 c ( 5) broadleaf shrubs FROM above type 6 (closed shrublands)
1061 c ( 6) dwarf trees (tundra) FROM above type 18 (tundra)
1062 c ( 7) bare soil FROM above type 16 (desert)
1063 c ( 8) light desert FROM above type 16 (desert)
1064 c ( 9) glacier FROM above type 15 (snow/ice)
1065 c ( 10) dark desert FROM above type 16 (desert)
1066 c (100) ocean FROM above type 17 (water)
1067 c
1068 c NOTE: snow-covered ground uses interpolated emissivities based on snow depth
1069 c =============================================================================
1070 c -----------------------------------------------------------------------------
1071 c Emmissivities for 12 bands in Fu/Liou
1072 c band 1: 4.5 - 5.3 um
1073 c band 2: 5.3 - 5.9 um
1074 c band 3: 5.9 - 7.1 um
1075 c band 4: 7.1 - 8.0 um
1076 c band 5: 8.0 - 9.1 um
1077 c band 6: 9.1 - 10.2 um
1078 c band 7: 10.2 - 12.5 um
1079 c band 8: 12.5 - 14.9 um
1080 c band 9: 14.9 - 18.5 um
1081 c band 10: 18.5 - 25.0 um
1082 c band 11: 25.0 - 35.7 um
1083 c band 12: 35.7 - oo um
1084 c
1085 c-------------------------------------------------------------------------
1086 data ((emis(i,j),i=1,12),j=1,11) /
1087 C evergreen needleleaf
1088 & 0.9891, 0.9892, 0.9900, 0.9914, 0.9908, 0.9903,
1089 & 0.9898, 0.9948, 1.0000, 1.0000, 1.0000, 1.0000,
1090 C deciduous needleleaf
1091 & 0.9849, 0.9856, 0.9841, 0.9831, 0.9789, 0.9805,
1092 & 0.9733, 0.9869, 1.0000, 1.0000, 1.0000, 1.0000,
1093 C evergreen needleleaf
1094 & 0.9891, 0.9892, 0.9900, 0.9914, 0.9908, 0.9903,
1095 & 0.9898, 0.9948, 1.0000, 1.0000, 1.0000, 1.0000,
1096 C grasslands
1097 & 0.9867, 0.9897, 0.9920, 0.9933, 0.9830, 0.9752,
1098 & 0.9853, 0.9928, 1.0000, 1.0000, 1.0000, 1.0000,
1099 C closed shrublands
1100 & 0.9490, 0.9697, 0.9738, 0.9712, 0.9474, 0.9582,
1101 & 0.9663, 0.9747, 0.9836, 0.9836, 0.9836, 0.9836,
1102 C tundra
1103 & 0.9469, 0.9670, 0.9883, 0.9795, 0.9751, 0.9767,
1104 & 0.9920, 0.9888, 0.9888, 0.9888, 0.9888, 0.9888,
1105 C barren
1106 & 0.8353, 0.9163, 0.9342, 0.9229, 0.8354, 0.8766,
1107 & 0.9210, 0.9262, 0.9345, 0.9345, 0.9345, 0.9345,
1108 C barren
1109 & 0.8353, 0.9163, 0.9342, 0.9229, 0.8354, 0.8766,
1110 & 0.9210, 0.9262, 0.9345, 0.9345, 0.9345, 0.9345,
1111 C snow/ice
1112 & 0.9998, 0.9998, 0.9998, 0.9998, 0.9998, 0.9999,
1113 & 0.9997, 0.9994, 0.9995, 0.9995, 0.9995, 0.9995,
1114 C barren
1115 & 0.8353, 0.9163, 0.9342, 0.9229, 0.8354, 0.8766,
1116 & 0.9210, 0.9262, 0.9345, 0.9345, 0.9345, 0.9345,
1117 C water
1118 & 0.9788, 0.9833, 0.9819, 0.9820, 0.9835, 0.9865,
1119 & 0.9886, 0.9719, 0.9719, 0.9719, 0.9719, 0.9719/
1120
1121 #include "snwmid.h"
1122
1123 c Convert to the 10 bands needed by Chou Radiation
1124 c ------------------------------------------------
1125 do i=1,numpts
1126
1127 c land points
1128 c------------
1129 if(ityp(i).le.10)then
1130 newemis(i, 1) = (emis( 1,ityp(i))+emis(2,ityp(i)))/2.
1131 newemis(i, 2) = (emis( 2,ityp(i))+emis(3,ityp(i)))/2.
1132 newemis(i, 3) = (emis( 4,ityp(i))+emis(5,ityp(i)))/2.
1133 newemis(i, 4) = emis( 6,ityp(i))
1134 newemis(i, 5) = emis( 7,ityp(i))
1135 newemis(i, 6) = emis( 8,ityp(i))
1136 newemis(i, 7) = emis( 9,ityp(i))
1137 newemis(i, 8) = (emis(10,ityp(i))+emis(11,ityp(i)))/2.
1138 newemis(i, 9) = emis(12,ityp(i))
1139 newemis(i,10) = emis( 4,ityp(i))
1140
1141 c modify emissivity for snow based on snow depth (like albedo)
1142 c-------------------------------------------------------------
1143 if(snowdepth (i).gt.0.) then
1144 fac = snowdepth(i) / (snowdepth(i) + snwmid(ityp(i)))
1145 newemis(i, 1) = newemis(i, 1) + (((emis( 1,9)+emis( 2,9))/2.)
1146 . - newemis(i, 1)) * fac
1147 newemis(i, 2) = newemis(i, 2) + (((emis( 2,9)+emis( 3,9))/2.)
1148 . - newemis(i, 2)) * fac
1149 newemis(i, 3) = newemis(i, 3) + (((emis( 4,9)+emis( 5,9))/2.)
1150 . - newemis(i, 3)) * fac
1151 newemis(i, 4) = newemis(i, 4) + (emis( 6,9)
1152 . - newemis(i, 4)) * fac
1153 newemis(i, 5) = newemis(i, 5) + (emis( 7,9)
1154 . - newemis(i, 5)) * fac
1155 newemis(i, 6) = newemis(i, 6) + (emis( 8,9)
1156 . - newemis(i, 6)) * fac
1157 newemis(i, 7) = newemis(i, 7) + (emis( 9,9)
1158 . - newemis(i, 7)) * fac
1159 newemis(i, 8) = newemis(i, 8) + (((emis(10,9)+emis(11,9))/2.)
1160 . - newemis(i, 8)) * fac
1161 newemis(i, 9) = newemis(i, 9) + (emis(12,9)
1162 . - newemis(i, 9)) * fac
1163 newemis(i,10) = newemis(i,10) + (emis( 4,9)
1164 . - newemis(i,10)) * fac
1165 endif
1166
1167 c open water
1168 c-----------
1169 else
1170 if(fraci(i).eq.0.)then
1171 newemis(i, 1) = (emis( 1,11)+emis(2,11))/2.
1172 newemis(i, 2) = (emis( 2,11)+emis(3,11))/2.
1173 newemis(i, 3) = (emis( 4,11)+emis(5,11))/2.
1174 newemis(i, 4) = emis( 6,11)
1175 newemis(i, 5) = emis( 7,11)
1176 newemis(i, 6) = emis( 8,11)
1177 newemis(i, 7) = emis( 9,11)
1178 newemis(i, 8) = (emis(10,11)+emis(11,11))/2.
1179 newemis(i, 9) = emis(12,11)
1180 newemis(i,10) = emis( 4,11)
1181
1182 c sea ice (like glacier and snow)
1183 c--------------------------------
1184 else
1185 newemis(i, 1) = (emis( 1,9)+emis(2,9))/2.
1186 newemis(i, 2) = (emis( 2,9)+emis(3,9))/2.
1187 newemis(i, 3) = (emis( 4,9)+emis(5,9))/2.
1188 newemis(i, 4) = emis( 6,9)
1189 newemis(i, 5) = emis( 7,9)
1190 newemis(i, 6) = emis( 8,9)
1191 newemis(i, 7) = emis( 9,9)
1192 newemis(i, 8) = (emis(10,9)+emis(11,9))/2.
1193 newemis(i, 9) = emis(12,9)
1194 newemis(i,10) = emis( 4,9)
1195 endif
1196 endif
1197 enddo
1198
1199 return
1200 end
1201 subroutine get_landfrac(im,jm,nSx,nSy,bi,bj,maxtyp,surftype,
1202 . tilefrac,frac)
1203 C***********************************************************************
1204 C Purpose
1205 C To compute the total fraction of land within a model grid-box
1206 C
1207 C***********************************************************************
1208 implicit none
1209
1210 integer im,jm,nSx,nSy,bi,bj,maxtyp
1211 integer surftype(im,jm,maxtyp,nSx,nSy)
1212 _RL tilefrac(im,jm,maxtyp,nSx,nSy)
1213 _RL frac(im,jm)
1214
1215 integer i,j,k
1216
1217 do j=1,jm
1218 do i=1,im
1219 frac(i,j) = 0.0
1220 enddo
1221 enddo
1222
1223 do k=1,maxtyp
1224 do j=1,jm
1225 do i=1,im
1226 if( (surftype(i,j,k,bi,bj).lt.100.).and.
1227 . (tilefrac(i,j,k,bi,bj).gt.0.0))then
1228 frac(i,j) = frac(i,j) + tilefrac(i,j,k,bi,bj)
1229 endif
1230 enddo
1231 enddo
1232 enddo
1233
1234 return
1235 end

  ViewVC Help
Powered by ViewVC 1.1.22