/[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.20 - (show annotations) (download)
Tue Jul 27 04:10:57 2004 UTC (19 years, 10 months ago) by molod
Branch: MAIN
CVS Tags: checkpoint54d_post
Changes since 1.19: +2 -4 lines
Debugging

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

  ViewVC Help
Powered by ViewVC 1.1.22