/[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.12 - (show annotations) (download)
Wed Jul 14 17:31:58 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.11: +3 -3 lines
Debugging

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

  ViewVC Help
Powered by ViewVC 1.1.22