/[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.15 - (show annotations) (download)
Fri Jul 23 22:32:28 2004 UTC (19 years, 11 months ago) by molod
Branch: MAIN
Changes since 1.14: +27 -20 lines
Debugging - got into main sequence!

1 C $Header: /u/gcmpack/MITgcm/pkg/fizhi/update_earth_exports.F,v 1.14 2004/07/16 20:08:08 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 "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 myTime, myIter, myThid
31
32 logical alarm
33 external alarm
34 real lats(sNx,sNy), lons(sNx,sNy), cosz(sNx,sNy)
35 real fraci(sNx,sNy), fracl(sNx,sNy)
36 real ficetile(nchp)
37 real radius
38 real tmpij(sNx,sNy)
39 real tmpchp(nchp)
40 integer i, j, n, bi, bj
41 integer im1, im2, jm1, jm2, idim1, idim2, jdim1, jdim2
42 integer sec, day, month
43 integer nmonf,ndayf,nsecf
44 nsecf(n) = n/10000*3600 + mod(n,10000)/100* 60 + mod(n,100)
45 nmonf(n) = mod(n,10000)/100
46 ndayf(n) = mod(n,100)
47
48 idim1 = 1-OLx
49 idim2 = sNx+OLx
50 jdim1 = 1-OLy
51 jdim2 = sNy+OLy
52 im1 = 1
53 im2 = sNx
54 jm1 = 1
55 jm2 = sNy
56 month = nmonf(nymd)
57 day = ndayf(nymd)
58 sec = nsecf(nhms)
59
60 do bj = myByLo(myThid), myByHi(myThid)
61 do bi = myBxLo(myThid), myBxHi(myThid)
62 do j = jm1,jm2
63 do i = im1,im2
64 lons(i,j) = xC(i,j,bi,bj)
65 lats(i,j) = yC(i,j,bi,bj)
66 enddo
67 enddo
68
69 call get_landfrac(im2,jm2,nSx,nSy,bi,bj,maxtyp,surftype,tilefrac,
70 . fracl)
71
72 do j = jm1,jm2
73 do i = im1,im2
74 if(sice(i,j,bi,bj).gt.0.) then
75 fraci(i,j) = 1.
76 else
77 fraci(i,j) = 0.
78 endif
79 enddo
80 enddo
81
82 C***********************************************************************
83 C* Get Leaf-Area-Index and Greenness Index *
84 C***********************************************************************
85
86 if( alarm('turb') .or. alarm('radsw') ) then
87 call getlgr (sec,month,day,chlt,ityp,nchpland,nchp,nSx,nSy,bi,bj,
88 . alai,agrn )
89 endif
90
91 C **********************************************************************
92 C Compute Surface Albedo
93 C **********************************************************************
94
95 if( alarm('radsw') ) then
96 call astro(nymd,nhms,lats,lons,im2*jm2,cosz,radius)
97 call getalb(sec,month,day,cosz,snodep,fraci,fracl,im2,jm2,nchp,
98 . nchptot,nchpland,nSx,nSy,bi,bj,igrd,ityp,chfr,chlt,alai,agrn,
99 . albvisdr,albvisdf,albnirdr,albnirdf )
100 endif
101
102
103 C **********************************************************************
104 C Compute Surface Emissivity
105 C **********************************************************************
106
107 if( alarm('radlw') ) then
108 call grd2msc(fraci,im2,jm2,igrd,ficetile,nchp,nchp)
109 call getemiss(fracl,im2,jm2,nchp,nSx,nSy,bi,bj,igrd,ityp,chfr,
110 . snodep,ficetile,emiss)
111 endif
112
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,nchp
125 tmpchp(i) = tcanopy(i,bi,bj)
126 enddo
127 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),tmpchp,
128 . nchp,nchptot,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 #include "CPP_EEOPTIONS.h"
165
166 INTEGER IRUN
167 real AVISDR (IRUN), ANIRDR (IRUN), AVISDF (IRUN), ANIRDF (IRUN)
168 _RL VLAI(IRUN),VGRN (IRUN), SNW(IRUN)
169 REAL ZTH(IRUN)
170 INTEGER ITYP (IRUN)
171
172 _RL ALVDRS, ALIDRS
173 _RL ALVDRDL, ALIDRDL
174 _RL ALVDRDD, ALIDRDD
175 _RL ALVDRI, ALIDRI
176 _RL minval
177 external minval
178
179 C Albedo of soil for visible direct solar radiation.
180 PARAMETER ( ALVDRS = 0.100 )
181 C Albedo of soil for infra-red direct solar radiation.
182 PARAMETER ( ALIDRS = 0.200 )
183 C Albedo of light desert for visible direct solar radiation.
184 PARAMETER ( ALVDRDL = 0.300 )
185 C Albedo of light desert for infra-red direct solar radiation.
186 PARAMETER ( ALIDRDL = 0.350 )
187 C Albedo of dark desert for visible direct solar radiation.
188 PARAMETER ( ALVDRDD = 0.250 )
189 C Albedo of dark desert for infra-red direct solar radiation.
190 PARAMETER ( ALIDRDD = 0.300 )
191 C Albedo of ice for visible direct solar radiation.
192 PARAMETER ( ALVDRI = 0.800 )
193 C Albedo of ice for infra-red direct solar radiation.
194 PARAMETER ( ALIDRI = 0.800 )
195
196 * --------------------------------------------------------------------------------------------
197
198 INTEGER NTYPS
199 INTEGER NLAI
200 _RL ZERO, ONE
201 _RL EPSLN, BLAI, DLAI
202 _RL ALATRM
203 PARAMETER (NLAI = 14 )
204 PARAMETER (EPSLN = 1.E-6)
205 PARAMETER (BLAI = 0.5)
206 PARAMETER (DLAI = 0.5)
207 PARAMETER (ZERO=0., ONE=1.0)
208 PARAMETER (ALATRM = BLAI + (NLAI - 1) * DLAI - EPSLN)
209 PARAMETER (NTYPS=10)
210
211
212 C ITYP: Vegetation type as follows:
213 C 1: BROADLEAF EVERGREEN TREES
214 C 2: BROADLEAF DECIDUOUS TREES
215 C 3: NEEDLELEAF TREES
216 C 4: GROUND COVER
217 C 5: BROADLEAF SHRUBS
218 C 6: DWARF TREES (TUNDRA)
219 C 7: BARE SOIL
220 C 8: LIGHT DESERT
221 C 9: GLACIER
222 C 10: DARK DESERT
223 C
224
225 INTEGER I, LAI
226 real FAC,GAMMA,BETA,ALPHA,DX,DY,ALA,GRN (2),SNWALB(4,NTYPS)
227 real COEFF
228
229 real ALVDR (NLAI, 2, NTYPS)
230 real BTVDR (NLAI, 2, NTYPS)
231 real GMVDR (NLAI, 2, NTYPS)
232 real ALIDR (NLAI, 2, NTYPS)
233 real BTIDR (NLAI, 2, NTYPS)
234 real 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 #ifdef CRAY
639 #ifdef 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 real coeff
684 real 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,nchpdim,
695 . nSx,nSy,bi,bj,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, nchpdim, nSx, nSy, bi, bj
707 _RL ALAI(nchpdim,nSx,nSy), AGRN(nchpdim,nSx,nSy)
708 _RL ALAT(nchpdim)
709 integer ITYP(nchpdim,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,nchptot,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,nchptot,nchpland,nSx,nSy,bi,bj
855 real 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 real alboc(im,jm)
872 real AVISDR(nchp),ANIRDR(nchp),AVISDF(nchp)
873 real ANIRDF(nchp)
874 real zenith(nchp)
875 real tmpij(im,jm)
876 integer i,j
877
878 DO I=1,IM
879 DO J=1,JM
880 ALBOC(I,J) = A0 + (A1 + (A2 + A3*cosz(I,J))*cosz(I,J))*cosz(I,J)
881 ALBVR(I,J,bi,bj) = ALBSI*FRACI(I,J) + ALBOC(I,J)*(ONE-FRACI(I,J))
882 ALBNR(I,J,bi,bj) = ALBVR(I,J,bi,bj)
883 ALBVF(I,J,bi,bj) = ALBSI * FRACI(I,J) + OCNALB * (ONE-FRACI(I,J))
884 ALBNF(I,J,bi,bj) = ALBVF(I,J,bi,bj)
885 ENDDO
886 ENDDO
887
888 C and now some conversions from grid space to tile space before sibalb
889
890 call grd2msc(cosz,im,jm,igrd,zenith,nchp,nchpland)
891
892 C and now call sibalb
893
894 call sibalb(avisdr,anirdr,avisdf,anirdf,alai(1,bi,bj),
895 . agrn(1,bi,bj),zenith,snodep(1,bi,bj),ityp(1,bi,bj),nchpland)
896
897 C finally some transformations back to grid space for albedos
898
899 print *,' In getalb, chfr: '
900 print *,(chfr(i,1,1),i=1,nchptot)
901
902 DO I=1,IM
903 DO J=1,JM
904 tmpij(i,j) = albvr(i,j,bi,bj)
905 ENDDO
906 ENDDO
907 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),avisdr,nchp,nchpland,
908 . fracg,tmpij,im,jm)
909
910 print *,' back from first msc2grd call '
911 stop
912
913 DO I=1,IM
914 DO J=1,JM
915 albvr(i,j,bi,bj) = tmpij(i,j)
916 ENDDO
917 ENDDO
918 DO I=1,IM
919 DO J=1,JM
920 tmpij(i,j) = albvf(i,j,bi,bj)
921 ENDDO
922 ENDDO
923 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),avisdf,nchp,nchpland,
924 . fracg,tmpij,im,jm)
925 DO I=1,IM
926 DO J=1,JM
927 albvf(i,j,bi,bj) = tmpij(i,j)
928 ENDDO
929 ENDDO
930 DO I=1,IM
931 DO J=1,JM
932 tmpij(i,j) = albnr(i,j,bi,bj)
933 ENDDO
934 ENDDO
935 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),anirdr,nchp,nchpland,
936 . fracg,tmpij,im,jm)
937 DO I=1,IM
938 DO J=1,JM
939 albnr(i,j,bi,bj) = tmpij(i,j)
940 ENDDO
941 ENDDO
942 DO I=1,IM
943 DO J=1,JM
944 tmpij(i,j) = albnf(i,j,bi,bj)
945 ENDDO
946 ENDDO
947 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),anirdf,nchp,nchpland,
948 . fracg,tmpij,im,jm)
949 DO I=1,IM
950 DO J=1,JM
951 albnf(i,j,bi,bj) = tmpij(i,j)
952 ENDDO
953 ENDDO
954
955 return
956 end
957
958 subroutine getemiss(fracg,im,jm,nchp,nSx,nSy,bi,bj,igrd,ityp,
959 . chfr,snowdep,fraci,emiss)
960 C***********************************************************************
961 C PURPOSE
962 C To act as an interface to routine to emissivity, which calculates
963 C ten bands of surface emissivities for use by the longwave radiation
964 C
965 C INPUT:
966 C fracg - real array in grid space of total land fraction [im,jm]
967 C im - model grid longitude dimension
968 C jm - model grid latitude dimension (number of lat. points)
969 C nchp - integer actual number of tiles in tile space
970 C nSx - number of processors in x-direction
971 C nSy - number of processors in y-direction
972 C bi - processors index in x-direction
973 C bj - processors index in y-direction
974 C igrd - integer array in tile space of grid point number for each
975 C tile [nchp]
976 C ityp - integer array in tile space of land surface type for each
977 C tile [nchp]
978 C chfr - real array in tile space of land surface type fraction for
979 C each tile [nchp]
980 C snowdep - real array in tile space of snow depth (liquid water equiv)
981 C in mm [nchp]
982 C fraci - real array in tile space of sea ice fraction [nchp]
983 C
984 C OUTPUT:
985 C emiss - real array [im,jm,10,nSx,nSy] - surface emissivity (frac)
986 C
987 C***********************************************************************
988 implicit none
989 #include "CPP_EEOPTIONS.h"
990 integer im,jm,nchp,nSx,nSy,bi,bj
991 real fracg(im,jm)
992 _RL chfr(nchp,nSx,nSy)
993 integer igrd(nchp,nSx,nSy), ityp(nchp,nSx,nSy)
994 _RL snowdep(nchp,nSx,nSy)
995 real fraci(nchp)
996 _RL emiss(im,jm,10,nSx,nSy)
997
998 real emisstile(nchp,10)
999 real tmpij(im,jm)
1000 integer i,j,k,n
1001
1002 do i = 1,10
1003 do n = 1,nchp
1004 emisstile(n,i) = 1.
1005 enddo
1006 enddo
1007
1008 c call emissivity to get values in tile space
1009 c -------------------------------------------
1010 call emissivity(snowdep(1,bi,bj),fraci,nchp,ityp(1,bi,bj),
1011 . emisstile)
1012
1013 c transform back to grid space for emissivities
1014 c ---------------------------------------------
1015 do k = 1,10
1016 do j = 1,jm
1017 do i = 1,im
1018 tmpij(i,j) = 0.0
1019 enddo
1020 enddo
1021 call msc2grd(igrd(1,bi,bj),chfr(1,bi,bj),emisstile(1,k),nchp,nchp,
1022 . fracg,tmpij,im,jm)
1023 do j = 1,jm
1024 do i = 1,im
1025 emiss(i,j,k,bi,bj) = tmpij(i,j)
1026 enddo
1027 enddo
1028 enddo
1029
1030 return
1031 end
1032
1033 subroutine emissivity (snowdepth,fraci,numpts,ityp,newemis)
1034 implicit none
1035 #include "CPP_EEOPTIONS.h"
1036 integer numpts
1037 integer ityp(numpts)
1038 _RL snowdepth(numpts)
1039 real fraci(numpts)
1040 real newemis(numpts,10)
1041
1042 real emis(12,11)
1043 real fac
1044 integer i,j
1045
1046 c-----------------------------------------------------------------------
1047 c NOTE: Emissivities were obtained for the following surface types:
1048 c ( 1) evergreen needleleaf = conifer
1049 c ( 2) evergreen broadleaf = conifer
1050 c ( 3) deciduous needleleaf = deciduous
1051 c ( 4) deciduous broadleaf = deciduous
1052 c ( 5) mixed forests = 1/2 conifer + 1/2 deciduous = tree
1053 c ( 6) closed shrublands = 3/4 tree + 1/4 quartz
1054 c ( 7) open shrubland = 1/4 tree + 3/4 quartz
1055 c ( 8) woody savannas = grass
1056 c ( 9) savannas = grass
1057 c (10) grasslands = grass
1058 c (11) permanent wetlands = 1/2 grass + 1/2 water
1059 c (12) croplands = grass
1060 c (13) urban = black body
1061 c (14) mosaic = 1/2 grass + 1/2 mixed forest
1062 c (15) snow/ice
1063 c (16) barren/sparsely vegetated = desert(quartz)
1064 c (17) water
1065 c (18) tundra = frost
1066 c
1067 c NOTE: Translation to Koster-Suarez surface types was as follows:
1068 c ( 1) broadleaf evergreen FROM above type 1 (conifer)
1069 c ( 2) broadleaf deciduous FROM above type 2 (deciduous)
1070 c ( 3) needleleaf evergreen FROM above type 1 (conifer)
1071 c ( 4) groundcover FROM above type 10 (grass)
1072 c ( 5) broadleaf shrubs FROM above type 6 (closed shrublands)
1073 c ( 6) dwarf trees (tundra) FROM above type 18 (tundra)
1074 c ( 7) bare soil FROM above type 16 (desert)
1075 c ( 8) light desert FROM above type 16 (desert)
1076 c ( 9) glacier FROM above type 15 (snow/ice)
1077 c ( 10) dark desert FROM above type 16 (desert)
1078 c (100) ocean FROM above type 17 (water)
1079 c
1080 c NOTE: snow-covered ground uses interpolated emissivities based on snow depth
1081 c =============================================================================
1082 c -----------------------------------------------------------------------------
1083 c Emmissivities for 12 bands in Fu/Liou
1084 c band 1: 4.5 - 5.3 um
1085 c band 2: 5.3 - 5.9 um
1086 c band 3: 5.9 - 7.1 um
1087 c band 4: 7.1 - 8.0 um
1088 c band 5: 8.0 - 9.1 um
1089 c band 6: 9.1 - 10.2 um
1090 c band 7: 10.2 - 12.5 um
1091 c band 8: 12.5 - 14.9 um
1092 c band 9: 14.9 - 18.5 um
1093 c band 10: 18.5 - 25.0 um
1094 c band 11: 25.0 - 35.7 um
1095 c band 12: 35.7 - oo um
1096 c
1097 c-------------------------------------------------------------------------
1098 data ((emis(i,j),i=1,12),j=1,11) /
1099 C evergreen needleleaf
1100 & 0.9891, 0.9892, 0.9900, 0.9914, 0.9908, 0.9903,
1101 & 0.9898, 0.9948, 1.0000, 1.0000, 1.0000, 1.0000,
1102 C deciduous needleleaf
1103 & 0.9849, 0.9856, 0.9841, 0.9831, 0.9789, 0.9805,
1104 & 0.9733, 0.9869, 1.0000, 1.0000, 1.0000, 1.0000,
1105 C evergreen needleleaf
1106 & 0.9891, 0.9892, 0.9900, 0.9914, 0.9908, 0.9903,
1107 & 0.9898, 0.9948, 1.0000, 1.0000, 1.0000, 1.0000,
1108 C grasslands
1109 & 0.9867, 0.9897, 0.9920, 0.9933, 0.9830, 0.9752,
1110 & 0.9853, 0.9928, 1.0000, 1.0000, 1.0000, 1.0000,
1111 C closed shrublands
1112 & 0.9490, 0.9697, 0.9738, 0.9712, 0.9474, 0.9582,
1113 & 0.9663, 0.9747, 0.9836, 0.9836, 0.9836, 0.9836,
1114 C tundra
1115 & 0.9469, 0.9670, 0.9883, 0.9795, 0.9751, 0.9767,
1116 & 0.9920, 0.9888, 0.9888, 0.9888, 0.9888, 0.9888,
1117 C barren
1118 & 0.8353, 0.9163, 0.9342, 0.9229, 0.8354, 0.8766,
1119 & 0.9210, 0.9262, 0.9345, 0.9345, 0.9345, 0.9345,
1120 C barren
1121 & 0.8353, 0.9163, 0.9342, 0.9229, 0.8354, 0.8766,
1122 & 0.9210, 0.9262, 0.9345, 0.9345, 0.9345, 0.9345,
1123 C snow/ice
1124 & 0.9998, 0.9998, 0.9998, 0.9998, 0.9998, 0.9999,
1125 & 0.9997, 0.9994, 0.9995, 0.9995, 0.9995, 0.9995,
1126 C barren
1127 & 0.8353, 0.9163, 0.9342, 0.9229, 0.8354, 0.8766,
1128 & 0.9210, 0.9262, 0.9345, 0.9345, 0.9345, 0.9345,
1129 C water
1130 & 0.9788, 0.9833, 0.9819, 0.9820, 0.9835, 0.9865,
1131 & 0.9886, 0.9719, 0.9719, 0.9719, 0.9719, 0.9719/
1132
1133 #include "snwmid.h"
1134
1135 c Convert to the 10 bands needed by Chou Radiation
1136 c ------------------------------------------------
1137 do i=1,numpts
1138
1139 c land points
1140 c------------
1141 if(ityp(i).le.10)then
1142 newemis(i, 1) = (emis( 1,ityp(i))+emis(2,ityp(i)))/2.
1143 newemis(i, 2) = (emis( 2,ityp(i))+emis(3,ityp(i)))/2.
1144 newemis(i, 3) = (emis( 4,ityp(i))+emis(5,ityp(i)))/2.
1145 newemis(i, 4) = emis( 6,ityp(i))
1146 newemis(i, 5) = emis( 7,ityp(i))
1147 newemis(i, 6) = emis( 8,ityp(i))
1148 newemis(i, 7) = emis( 9,ityp(i))
1149 newemis(i, 8) = (emis(10,ityp(i))+emis(11,ityp(i)))/2.
1150 newemis(i, 9) = emis(12,ityp(i))
1151 newemis(i,10) = emis( 4,ityp(i))
1152
1153 c modify emissivity for snow based on snow depth (like albedo)
1154 c-------------------------------------------------------------
1155 if(snowdepth (i).gt.0.) then
1156 fac = snowdepth(i) / (snowdepth(i) + snwmid(ityp(i)))
1157 newemis(i, 1) = newemis(i, 1) + (((emis( 1,9)+emis( 2,9))/2.)
1158 . - newemis(i, 1)) * fac
1159 newemis(i, 2) = newemis(i, 2) + (((emis( 2,9)+emis( 3,9))/2.)
1160 . - newemis(i, 2)) * fac
1161 newemis(i, 3) = newemis(i, 3) + (((emis( 4,9)+emis( 5,9))/2.)
1162 . - newemis(i, 3)) * fac
1163 newemis(i, 4) = newemis(i, 4) + (emis( 6,9)
1164 . - newemis(i, 4)) * fac
1165 newemis(i, 5) = newemis(i, 5) + (emis( 7,9)
1166 . - newemis(i, 5)) * fac
1167 newemis(i, 6) = newemis(i, 6) + (emis( 8,9)
1168 . - newemis(i, 6)) * fac
1169 newemis(i, 7) = newemis(i, 7) + (emis( 9,9)
1170 . - newemis(i, 7)) * fac
1171 newemis(i, 8) = newemis(i, 8) + (((emis(10,9)+emis(11,9))/2.)
1172 . - newemis(i, 8)) * fac
1173 newemis(i, 9) = newemis(i, 9) + (emis(12,9)
1174 . - newemis(i, 9)) * fac
1175 newemis(i,10) = newemis(i,10) + (emis( 4,9)
1176 . - newemis(i,10)) * fac
1177 endif
1178
1179 c open water
1180 c-----------
1181 else
1182 if(fraci(i).eq.0.)then
1183 newemis(i, 1) = (emis( 1,11)+emis(2,11))/2.
1184 newemis(i, 2) = (emis( 2,11)+emis(3,11))/2.
1185 newemis(i, 3) = (emis( 4,11)+emis(5,11))/2.
1186 newemis(i, 4) = emis( 6,11)
1187 newemis(i, 5) = emis( 7,11)
1188 newemis(i, 6) = emis( 8,11)
1189 newemis(i, 7) = emis( 9,11)
1190 newemis(i, 8) = (emis(10,11)+emis(11,11))/2.
1191 newemis(i, 9) = emis(12,11)
1192 newemis(i,10) = emis( 4,11)
1193
1194 c sea ice (like glacier and snow)
1195 c--------------------------------
1196 else
1197 newemis(i, 1) = (emis( 1,9)+emis(2,9))/2.
1198 newemis(i, 2) = (emis( 2,9)+emis(3,9))/2.
1199 newemis(i, 3) = (emis( 4,9)+emis(5,9))/2.
1200 newemis(i, 4) = emis( 6,9)
1201 newemis(i, 5) = emis( 7,9)
1202 newemis(i, 6) = emis( 8,9)
1203 newemis(i, 7) = emis( 9,9)
1204 newemis(i, 8) = (emis(10,9)+emis(11,9))/2.
1205 newemis(i, 9) = emis(12,9)
1206 newemis(i,10) = emis( 4,9)
1207 endif
1208 endif
1209 enddo
1210
1211 return
1212 end
1213 subroutine get_landfrac(im,jm,nSx,nSy,bi,bj,maxtyp,surftype,
1214 . tilefrac,frac)
1215 C***********************************************************************
1216 C Purpose
1217 C To compute the total fraction of land within a model grid-box
1218 C
1219 C***********************************************************************
1220 implicit none
1221 #include "CPP_EEOPTIONS.h"
1222
1223 integer im,jm,nSx,nSy,bi,bj,maxtyp
1224 integer surftype(im,jm,maxtyp,nSx,nSy)
1225 _RL tilefrac(im,jm,maxtyp,nSx,nSy)
1226 real frac(im,jm)
1227
1228 integer i,j,k
1229
1230 do j=1,jm
1231 do i=1,im
1232 frac(i,j) = 0.0
1233 enddo
1234 enddo
1235
1236 do k=1,maxtyp
1237 do j=1,jm
1238 do i=1,im
1239 if( (surftype(i,j,k,bi,bj).lt.100.).and.
1240 . (tilefrac(i,j,k,bi,bj).gt.0.0))then
1241 frac(i,j) = frac(i,j) + tilefrac(i,j,k,bi,bj)
1242 endif
1243 enddo
1244 enddo
1245 enddo
1246
1247 return
1248 end

  ViewVC Help
Powered by ViewVC 1.1.22