/[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.10 - (show annotations) (download)
Wed Jun 16 19:19:49 2004 UTC (20 years ago) by molod
Branch: MAIN
Changes since 1.9: +26 -18 lines
Developing

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

  ViewVC Help
Powered by ViewVC 1.1.22