/[MITgcm]/MITgcm/pkg/aim_v23/phy_radiat.F
ViewVC logotype

Contents of /MITgcm/pkg/aim_v23/phy_radiat.F

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


Revision 1.2 - (show annotations) (download)
Thu Mar 11 14:33:19 2004 UTC (20 years, 3 months ago) by jmc
Branch: MAIN
CVS Tags: checkpoint53, checkpoint53d_post, checkpoint52m_post, checkpoint53c_post, checkpoint53a_post, checkpoint52l_post, checkpoint52n_post, checkpoint53b_pre, checkpoint53b_post, checkpoint53d_pre
Changes since 1.1: +12 -6 lines
a) Treat separately land / ocean / sea-ice surface fluxes
   to allow implicit computation of land & sea-ice surface temp.
b) add snow precipitation
c) other (little) modifications for new land formulation.

1 C $Header: /u/gcmpack/MITgcm/pkg/aim_v23/phy_radiat.F,v 1.1 2002/11/22 17:17:03 jmc Exp $
2 C $Name: $
3
4 #include "AIM_OPTIONS.h"
5
6 SUBROUTINE SOL_OZ (SOLC, TYEAR, SLAT, CLAT,
7 O FSOL, OZONE, OZUPP, ZENIT, STRATZ,
8 I bi, bj, myThid)
9
10 C--
11 C-- SUBROUTINE SOL_OZ (SOLC,TYEAR)
12 C--
13 C-- Purpose: Compute the flux of incoming solar radiation
14 C-- and a climatological ozone profile for SW absorption
15 C-- Input: SOLC = solar constant (area averaged)
16 C-- TYEAR = time as fraction of year (0-1, 0 = 1jan.h00)
17 C SLAT = sin(lat)
18 C CLAT = cos(lat)
19 C-- Output: FSOL = flux of incoming solar radiation
20 C-- OZONE = flux absorbed by ozone (lower stratos.)
21 C-- OZUPP = flux absorbed by ozone (upper stratos.)
22 C-- ZENIT = function of solar zenith angle
23 C-- STRATZ = ?
24 C-- Updated common blocks: RADZON
25 C--
26
27 IMPLICIT NONE
28
29 C Resolution parameters
30
31 C-- size for MITgcm & Physics package :
32 #include "AIM_SIZE.h"
33
34 #include "EEPARAMS.h"
35
36 C Constants + functions of sigma and latitude
37 #include "com_physcon.h"
38
39 C Radiation constants
40 #include "com_radcon.h"
41
42 C-- Routine arguments:
43 INTEGER bi, bj, myThid
44 _RL SOLC, TYEAR
45 _RL SLAT(NGP), CLAT(NGP)
46 _RL FSOL(NGP), OZONE(NGP), OZUPP(NGP), ZENIT(NGP), STRATZ(NGP)
47
48 #ifdef ALLOW_AIM
49
50 C-- Local variables:
51 INTEGER J, NZEN
52
53 C- jmc: declare all local variables:
54 _RL ALPHA, CSR1, CSR2, COZ1, COZ2
55 _RL AZEN, RZEN, CZEN, SZEN, AST, FS0, FLAT2
56 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
57
58 C ALPHA = year phase ( 0 - 2pi, 0 = winter solstice = 22dec.h00 )
59 ALPHA = 4. _d 0*ASIN(1. _d 0)*(TYEAR+10. _d 0/365. _d 0)
60
61 CSR1=-0.796 _d 0*COS(ALPHA)
62 CSR2= 0.147 _d 0*COS(2. _d 0*ALPHA)-0.477 _d 0
63 COZ1= 1.0 _d 0 * COS(ALPHA)
64 COZ2= 1.8 _d 0
65 C
66 AZEN=1.0
67 NZEN=2
68
69 RZEN=-COS(ALPHA)*23.45 _d 0*ASIN(1. _d 0)/90. _d 0
70 CZEN=COS(RZEN)
71 SZEN=SIN(RZEN)
72
73 AST=0.025 _d 0
74 FS0=10. _d 0
75 C FS0=16.-8.*COS(ALPHA)
76
77 DO J=1,NGP
78
79 FLAT2 = 1.5 _d 0*SLAT(J)**2 - 0.5 _d 0
80
81 C solar radiation at the top
82 FSOL(J) = SOLC*
83 & MAX( 0. _d 0, 1. _d 0+CSR1*SLAT(J)+CSR2*FLAT2 )
84
85 C ozone depth in upper and lower stratosphere
86 OZUPP(J) = EPSSW*(1.-FLAT2)
87 OZONE(J) = EPSSW*(1.+COZ1*SLAT(J)+COZ2*FLAT2)
88
89 C zenith angle correction to (downward) absorptivity
90 ZENIT(J) = 1. + AZEN*
91 & (1. _d 0-(CLAT(J)*CZEN+SLAT(J)*SZEN))**NZEN
92
93 C ozone absorption in upper and lower stratosphere
94 OZUPP(J)=FSOL(J)*OZUPP(J)*ZENIT(J)
95 OZONE(J)=FSOL(J)*OZONE(J)*ZENIT(J)
96 STRATZ(J)=AST*FSOL(J)*CLAT(J)**3
97 & +MAX( FS0-FSOL(J), 0. _d 0 )
98
99 ENDDO
100
101 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
102
103 #endif /* ALLOW_AIM */
104 RETURN
105 END
106
107
108 SUBROUTINE RADSW (PSA,dpFac,QA,RH,ALB,
109 I FSOL, OZONE, OZUPP, ZENIT, STRATZ,
110 O TAU2, STRATC,
111 O ICLTOP,CLOUDC,FTOP,FSFC,DFABS,
112 I kGrd,bi,bj,myThid)
113 C--
114 C-- SUBROUTINE RADSW (PSA,QA,RH,ALB,
115 C-- & ICLTOP,CLOUDC,FTOP,FSFC,DFABS)
116 C--
117 C-- Purpose: Compute the absorption of shortwave radiation and
118 C-- initialize arrays for longwave-radiation routines
119 C-- Input: PSA = norm. surface pressure [p/p0] (2-dim)
120 C dpFac = cell delta_P fraction (3-dim)
121 C-- QA = specific humidity [g/kg] (3-dim)
122 C-- RH = relative humidity (3-dim)
123 C-- ALB = surface albedo (2-dim)
124 C-- Output: ICLTOP = cloud top level (2-dim)
125 C-- CLOUDC = total cloud cover (2-dim)
126 C-- FTOP = net downw. flux of sw rad. at the atm. top (2-dim)
127 C-- FSFC = net downw. flux of sw rad. at the surface (2-dim)
128 C-- DFABS = flux of sw rad. absorbed by each atm. layer (3-dim)
129 C Input: kGrd = Ground level index (2-dim)
130 C bi,bj = tile index
131 C myThid = Thread number for this instance of the routine
132 C--
133
134 IMPLICIT NONE
135
136 C Resolution parameters
137
138 C-- size for MITgcm & Physics package :
139 #include "AIM_SIZE.h"
140
141 #include "EEPARAMS.h"
142
143 C Constants + functions of sigma and latitude
144
145 #include "com_physcon.h"
146
147 C Radiation parameters
148
149 #include "com_radcon.h"
150
151 C-- Routine arguments:
152 _RL PSA(NGP),dpFac(NGP,NLEV),QA(NGP,NLEV),RH(NGP,NLEV)
153 _RL ALB(NGP,0:3)
154 INTEGER ICLTOP(NGP)
155 _RL CLOUDC(NGP), FTOP(NGP), FSFC(NGP,0:3), DFABS(NGP,NLEV)
156
157 _RL FSOL(NGP), OZONE(NGP), OZUPP(NGP), ZENIT(NGP), STRATZ(NGP)
158 _RL TAU2(NGP,NLEV,NBAND),STRATC(NGP)
159 c _RL FLUX(NGP,4)
160
161 INTEGER kGrd(NGP)
162 INTEGER bi, bj, myThid
163
164 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
165
166 #ifdef ALLOW_AIM
167
168 C-- Local variables:
169 c _RL QCLOUD(NGP), ACLOUD(NGP), PSAZ(NGP),
170 _RL QCLOUD(NGP), ACLOUD(NGP),
171 & ALBTOP(NGP,NLEV), FREFL(NGP,NLEV), FLUX(NGP,2)
172
173 C- jmc: define "FLUX" as a local variable & remove Equivalences:
174 c EQUIVALENCE (ALBTOP(1,1),TAU2(1,1,3))
175 c EQUIVALENCE ( FREFL(1,1),TAU2(1,1,4))
176
177 INTEGER NL1(NGP)
178 INTEGER K, J
179
180 C- jmc: declare local variables:
181 _RL FBAND1, FBAND2, RRCL, RQCL, DQACL, QACL3
182 _RL ABS1, DELTAP
183 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
184
185 FBAND2=0.05 _d 0
186 FBAND1=1.-FBAND2
187
188 DO J=1,NGP
189 NL1(J)=kGrd(J)-1
190 ENDDO
191
192 C-- 1. Cloud cover:
193 C defined as a linear fun. of the maximum relative humidity
194 C in all tropospheric layers above PBL:
195 C CLOUDC = 0 for RHmax < RHCL1, = 1 for RHmax > RHCL2.
196 C This value is reduced by a factor (Qbase/QACL) if the
197 C cloud-base absolute humidity Qbase < QACL.
198 C
199
200 RRCL=1./(RHCL2-RHCL1)
201 RQCL=1./QACL2
202 C
203 DO J=1,NGP
204 CLOUDC(J)=0.
205 QCLOUD(J)=0.
206 IF (kGrd(J).NE.0)
207 & QCLOUD(J)= MAX( QA(J,kGrd(J)), QA(J,NL1(J)) )
208 ICLTOP(J)= kGrd(J)
209 FREFL(J,1)=0.
210 ENDDO
211
212 DO K=1,NLEV
213 DO J=1,NGP
214 ALBTOP(J,K)=0.
215 ENDDO
216 ENDDO
217 C
218 DQACL=(QACL2-QACL1)/(0.5 _d 0 - SIG(2))
219 DO J=1,NGP
220 DO K=NL1(J),2,-1
221 QACL3=MIN(QACL2,QACL1+DQACL*(SIG(K)-SIG(2)))
222 IF (RH(J,K).GT.RHCL1.AND.QA(J,K).GT.QACL1) THEN
223 CLOUDC(J)=MAX(CLOUDC(J),RH(J,K)-RHCL1)
224 IF (QA(J,K).GT.QACL3) ICLTOP(J)=K
225 ENDIF
226 ENDDO
227 ENDDO
228
229 DO J=1,NGP
230 CLOUDC(J)=MIN(1. _d 0,CLOUDC(J)*RRCL)
231 IF (CLOUDC(J).GT.0.0) THEN
232 CLOUDC(J)=CLOUDC(J)*MIN(1. _d 0,QCLOUD(J)*RQCL)
233 ALBTOP(J,ICLTOP(J))=ALBCL*CLOUDC(J)
234 ELSE
235 ICLTOP(J)=NLEV+1
236 ENDIF
237 ENDDO
238
239 C
240 C-- 2. Shortwave transmissivity:
241 C function of layer mass, ozone (in the statosphere),
242 C abs. humidity and cloud cover (in the troposphere)
243
244 DO J=1,NGP
245 c_FM PSAZ(J)=PSA(J)*ZENIT(J)
246 ACLOUD(J)=CLOUDC(J)*(ABSCL1+ABSCL2*QCLOUD(J))
247 ENDDO
248
249 DO J=1,NGP
250 c_FM DELTAP=PSAZ(J)*DSIG(1)
251 DELTAP=ZENIT(J)*DSIG(1)*dpFac(J,1)
252 TAU2(J,1,1)=EXP(-DELTAP*ABSDRY)
253 ENDDO
254 C
255 DO J=1,NGP
256 DO K=2,NL1(J)
257 c_FM ABS1=ABSDRY+ABSAER*SIG(K)**2
258 c_FM DELTAP=PSAZ(J)*DSIG(K)
259 ABS1=ABSDRY+ABSAER*(SIG(K)/PSA(J))**2
260 DELTAP=ZENIT(J)*DSIG(K)*dpFac(J,K)
261 IF (K.EQ.ICLTOP(J)) THEN
262 TAU2(J,K,1)=EXP(-DELTAP*
263 & (ABS1+ABSWV1*QA(J,K)+2.*ACLOUD(J)))
264 ELSE IF (K.GT.ICLTOP(J)) THEN
265 TAU2(J,K,1)=EXP(-DELTAP*
266 & (ABS1+ABSWV1*QA(J,K)+ACLOUD(J)))
267 ELSE
268 TAU2(J,K,1)=EXP(-DELTAP*(ABS1+ABSWV1*QA(J,K)))
269 ENDIF
270 ENDDO
271 ENDDO
272
273 c_FM ABS1=ABSDRY+ABSAER*SIG(NLEV)**2
274 DO J=1,NGP
275 K = kGrd(J)
276 ABS1=ABSDRY+ABSAER*(SIG(K)/PSA(J))**2
277 c_FM DELTAP=PSAZ(J)*DSIG(NLEV)
278 DELTAP=ZENIT(J)*DSIG(K)*dpFac(J,K)
279 TAU2(J,K,1)=EXP(-DELTAP*(ABS1+ABSWV1*QA(J,K)))
280 ENDDO
281
282 DO J=1,NGP
283 DO K=2,kGrd(J)
284 DELTAP=ZENIT(J)*DSIG(K)*dpFac(J,K)
285 TAU2(J,K,2)=EXP(-DELTAP*ABSWV2*QA(J,K))
286 ENDDO
287 ENDDO
288 C
289 C--- 3. Shortwave downward flux
290 C
291 C 3.1 Absorption in the stratosphere
292
293 C 3.1.1 Initialization of fluxes (subtracting
294 C ozone absorption in the upper stratosphere)
295
296 DO J=1,NGP
297 FTOP(J) =FSOL(J)
298 FLUX(J,1)=FSOL(J)*FBAND1-OZUPP(J)
299 FLUX(J,2)=FSOL(J)*FBAND2
300 STRATC(J)=STRATZ(J)*PSA(J)
301 ENDDO
302
303 C 3.1.2 Ozone and dry-air absorption
304 C in the lower (modelled) stratosphere
305
306 DO J=1,NGP
307 DFABS(J,1)=FLUX(J,1)
308 FLUX (J,1)=TAU2(J,1,1)*(FLUX(J,1)-OZONE(J)*PSA(J))
309 DFABS(J,1)=DFABS(J,1)-FLUX(J,1)
310 ENDDO
311
312 C 3.3 Absorption and reflection in the troposphere
313 C
314 DO J=1,NGP
315 DO K=2,kGrd(J)
316 FREFL(J,K)=FLUX(J,1)*ALBTOP(J,K)
317 FLUX (J,1)=FLUX(J,1)-FREFL(J,K)
318 DFABS(J,K)=FLUX(J,1)
319 FLUX (J,1)=TAU2(J,K,1)*FLUX(J,1)
320 DFABS(J,K)=DFABS(J,K)-FLUX(J,1)
321 ENDDO
322 ENDDO
323
324 DO J=1,NGP
325 DO K=2,kGrd(J)
326 DFABS(J,K)=DFABS(J,K)+FLUX(J,2)
327 FLUX (J,2)=TAU2(J,K,2)*FLUX(J,2)
328 DFABS(J,K)=DFABS(J,K)-FLUX(J,2)
329 ENDDO
330 ENDDO
331
332 C
333 C--- 4. Shortwave upward flux
334 C
335 C 4.1 Absorption and reflection at the surface
336 C
337 DO J=1,NGP
338 C for each surface type:
339 FSFC(J,1)=FLUX(J,1)*(1.-ALB(J,1))+FLUX(J,2)
340 FSFC(J,2)=FLUX(J,1)*(1.-ALB(J,2))+FLUX(J,2)
341 FSFC(J,3)=FLUX(J,1)*(1.-ALB(J,3))+FLUX(J,2)
342 C weighted average according to land/sea/sea-ice fraction:
343 FSFC(J,0)=FLUX(J,1)+FLUX(J,2)
344 FLUX(J,1)=FLUX(J,1)*ALB(J,0)
345 FSFC(J,0)=FSFC(J,0)-FLUX(J,1)
346 ENDDO
347 C
348 C 4.2 Absorption of upward flux
349 C
350 DO J=1,NGP
351 DO K=kGrd(J),1,-1
352 DFABS(J,K)=DFABS(J,K)+FLUX(J,1)
353 FLUX (J,1)=TAU2(J,K,1)*FLUX(J,1)
354 DFABS(J,K)=DFABS(J,K)-FLUX(J,1)
355 FLUX (J,1)=FLUX(J,1)+FREFL(J,K)
356 ENDDO
357 ENDDO
358 C
359 C 4.3 Net solar radiation = incoming - outgoing
360 C
361 DO J=1,NGP
362 FTOP(J)=FTOP(J)-FLUX(J,1)
363 ENDDO
364
365 C
366 C--- 5. Initialization of longwave radiation model
367 C
368 C 5.1 Longwave transmissivity:
369 C function of layer mass, abs. humidity and cloud cover.
370
371 DO J=1,NGP
372 ACLOUD(J)=CLOUDC(J)*(ABLCL1+ABLCL2*QCLOUD(J))
373 ENDDO
374
375 DO J=1,NGP
376 c_FM DELTAP=PSA(J)*DSIG(1)
377 DELTAP=DSIG(1)*dpFac(J,1)
378 TAU2(J,1,1)=EXP(-DELTAP*ABLWIN)
379 TAU2(J,1,2)=EXP(-DELTAP*ABLCO2)
380 TAU2(J,1,3)=1.
381 TAU2(J,1,4)=1.
382 ENDDO
383
384 DO K=2,NLEV
385 DO J=1,NGP
386 c_FM DELTAP=PSA(J)*DSIG(K)
387 DELTAP=DSIG(K)*dpFac(J,K)
388 IF ( K.GE.ICLTOP(J).AND.K.NE.kGrd(J) ) THEN
389 TAU2(J,K,1)=EXP(-DELTAP*(ABLWIN+ACLOUD(J)))
390 ELSE
391 TAU2(J,K,1)=EXP(-DELTAP*ABLWIN)
392 ENDIF
393 TAU2(J,K,2)=EXP(-DELTAP*ABLCO2)
394 TAU2(J,K,3)=EXP(-DELTAP*ABLWV1*QA(J,K))
395 TAU2(J,K,4)=EXP(-DELTAP*ABLWV2*QA(J,K))
396 ENDDO
397 ENDDO
398 C
399 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
400 #endif /* ALLOW_AIM */
401
402 RETURN
403 END
404
405
406 SUBROUTINE RADLW (IMODE,TA,TS,ST4S,
407 I OZUPP, STRATC, TAU2,
408 & FLUX, ST4A,
409 & FTOP,FSFC,DFABS,
410 I kGrd,bi,bj,myThid)
411 C--
412 C-- SUBROUTINE RADLW (IMODE,TA,TS,ST4S,
413 C-- & FTOP,FSFC,DFABS)
414 C--
415 C-- Purpose: Compute the absorption of longwave radiation
416 C-- Input: IMODE = index for operation mode
417 C-- -1 : downward flux only
418 C-- 0 : downward + upward flux
419 C-- +1 : upward flux only
420 C-- TA = absolute temperature (3-dim)
421 C-- TS = surface temperature [if IMODE=0,1]
422 C-- ST4S = surface blackbody emission [if IMODE=1]
423 C-- FSFC = FSFC output from RADLW(-1,... ) [if IMODE=1]
424 C-- DFABS = DFABS output from RADLW(-1,... ) [if IMODE=1]
425 C-- Output: ST4S = surface blackbody emission [if IMODE=0]
426 C-- FTOP = outgoing flux of lw rad. at the top [if IMODE=0,1]
427 C-- FSFC = downward flux of lw rad. at the sfc. [if IMODE= -1]
428 C-- net upw. flux of lw rad. at the sfc. [if IMODE=0,1]
429 C-- DFABS = flux of lw rad. absorbed by each atm. layer (3-dim)
430 C Input: kGrd = Ground level index (2-dim)
431 C bi,bj = tile index
432 C myThid = Thread number for this instance of the routine
433 C--
434
435 IMPLICIT NONE
436
437 C Resolution parameters
438
439 C-- size for MITgcm & Physics package :
440 #include "AIM_SIZE.h"
441
442 #include "EEPARAMS.h"
443
444 C Number of radiation bands with tau < 1
445 c INTEGER NBAND
446 c PARAMETER ( NBAND=4 )
447
448 C Constants + functions of sigma and latitude
449 #include "com_physcon.h"
450
451 C Radiation parameters
452 #include "com_radcon.h"
453
454 C-- Routine arguments:
455 INTEGER IMODE
456 _RL TA(NGP,NLEV), TS(NGP), ST4S(NGP)
457 _RL FTOP(NGP), FSFC(NGP), DFABS(NGP,NLEV)
458 _RL OZUPP(NGP), STRATC(NGP)
459 _RL TAU2(NGP,NLEV,NBAND), FLUX(NGP,NBAND), ST4A(NGP,NLEV,2)
460
461 INTEGER kGrd(NGP)
462 INTEGER bi,bj,myThid
463
464 #ifdef ALLOW_AIM
465
466 C-- Local variables:
467 INTEGER K, J, JB
468 c INTEGER J0, Jl, I2
469 INTEGER NL1(NGP)
470
471 C- jmc: declare all local variables:
472 _RL REFSFC, BRAD, EMIS
473 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
474
475 DO J=1,NGP
476 NL1(J)=kGrd(J)-1
477 ENDDO
478
479 REFSFC=1.-EMISFC
480
481 IF (IMODE.EQ.1) GO TO 410
482
483 C--- 1. Blackbody emission from atmospheric full and half levels.
484 C Temperature is interpolated as a linear function of ln sigma.
485 C At the lower boundary, the emission is linearly extrapolated;
486 C at the upper boundary, the atmosphere is assumed isothermal.
487
488 DO K=1,NLEV
489 DO J=1,NGP
490 ST4A(J,K,1)=TA(J,K)*TA(J,K)
491 ST4A(J,K,1)=SBC*ST4A(J,K,1)*ST4A(J,K,1)
492 ENDDO
493 ENDDO
494 C
495 DO K=1,NLEV-1
496 DO J=1,NGP
497 ST4A(J,K,2)=TA(J,K)+WVI(K,2)*(TA(J,K+1)-TA(J,K))
498 ST4A(J,K,2)=ST4A(J,K,2)*ST4A(J,K,2)
499 ST4A(J,K,2)=SBC*ST4A(J,K,2)*ST4A(J,K,2)
500 ENDDO
501 ENDDO
502 C
503 DO J=1,NGP
504 c ST4A(J,NLEV,2)=ST4A(J,NLEV,1)
505 K=kGrd(J)
506 ST4A(J,K,2)=2.*ST4A(J,K,1)-ST4A(J,NL1(J),2)
507 ENDDO
508
509 C--- 2. Initialization
510 C--- (including the stratospheric correction term)
511
512 DO J=1,NGP
513 FTOP(J) = 0.
514 FSFC(J) = STRATC(J)
515 DFABS(J,1)=-STRATC(J)
516 ENDDO
517
518 DO K=2,NLEV
519 DO J=1,NGP
520 DFABS(J,K)=0.
521 ENDDO
522 ENDDO
523
524 C--- 3. Emission ad absorption of longwave downward flux.
525 C Downward emission is an average of the emission from the full level
526 C and the half-level below, weighted according to the transmissivity
527 C of the layer.
528
529 C 3.1 Stratosphere
530
531 K=1
532 DO JB=1,2
533 DO J=1,NGP
534 BRAD=ST4A(J,K,2)+TAU2(J,K,JB)*(ST4A(J,K,1)-ST4A(J,K,2))
535 EMIS=FBAND(NINT(TA(J,K)),JB)*(1.-TAU2(J,K,JB))
536 FLUX(J,JB)=EMIS*BRAD
537 DFABS(J,K)=DFABS(J,K)-FLUX(J,JB)
538 ENDDO
539 ENDDO
540
541 DO JB=3,NBAND
542 DO J=1,NGP
543 FLUX(J,JB)=0.
544 ENDDO
545 ENDDO
546
547 C 3.2 Troposphere
548
549 DO JB=1,NBAND
550 DO J=1,NGP
551 DO K=2,kGrd(J)
552 BRAD=ST4A(J,K,2)+TAU2(J,K,JB)*(ST4A(J,K,1)-ST4A(J,K,2))
553 EMIS=FBAND(NINT(TA(J,K)),JB)*(1.-TAU2(J,K,JB))
554 DFABS(J,K)=DFABS(J,K)+FLUX(J,JB)
555 FLUX(J,JB)=TAU2(J,K,JB)*FLUX(J,JB)+EMIS*BRAD
556 DFABS(J,K)=DFABS(J,K)-FLUX(J,JB)
557 ENDDO
558 ENDDO
559 ENDDO
560
561 DO JB=1,NBAND
562 DO J=1,NGP
563 FSFC(J)=FSFC(J)+EMISFC*FLUX(J,JB)
564 ENDDO
565 ENDDO
566
567 IF (IMODE.EQ.-1) RETURN
568
569 C--- 4. Emission ad absorption of longwave upward flux
570 C Upward emission is an average of the emission from the full level
571 C and the half-level above, weighted according to the transmissivity
572 C of the layer (for the top layer, full-level emission is used).
573 C Surface lw emission in "band 0" goes directly into FTOP.
574
575 C 4.1 Surface
576
577 DO J=1,NGP
578 ST4S(J)=TS(J)*TS(J)
579 ST4S(J)=SBC*ST4S(J)*ST4S(J)
580 ENDDO
581
582 C Entry point for upward-only mode (IMODE=1)
583 410 CONTINUE
584
585 DO J=1,NGP
586 ST4S(J)=EMISFC*ST4S(J)
587 FSFC(J)=ST4S(J)-FSFC(J)
588 FTOP(J)=FTOP(J)+FBAND(NINT(TS(J)),0)*ST4S(J)
589 ENDDO
590
591 DO JB=1,NBAND
592 DO J=1,NGP
593 FLUX(J,JB)=FBAND(NINT(TS(J)),JB)*ST4S(J)
594 & +REFSFC*FLUX(J,JB)
595 ENDDO
596 ENDDO
597
598 C 4.2 Troposphere
599
600 DO JB=1,NBAND
601 DO J=1,NGP
602 DO K=kGrd(J),2,-1
603 BRAD=ST4A(J,K-1,2)+TAU2(J,K,JB)*(ST4A(J,K,1)-ST4A(J,K-1,2))
604 EMIS=FBAND(NINT(TA(J,K)),JB)*(1.-TAU2(J,K,JB))
605 DFABS(J,K)=DFABS(J,K)+FLUX(J,JB)
606 FLUX(J,JB)=TAU2(J,K,JB)*FLUX(J,JB)+EMIS*BRAD
607 DFABS(J,K)=DFABS(J,K)-FLUX(J,JB)
608 ENDDO
609 ENDDO
610 ENDDO
611
612 C 4.3 Stratosphere
613
614 K=1
615 DO JB=1,2
616 DO J=1,NGP
617 EMIS=FBAND(NINT(TA(J,K)),JB)*(1.-TAU2(J,K,JB))
618 DFABS(J,K)=DFABS(J,K)+FLUX(J,JB)
619 FLUX(J,JB)=TAU2(J,K,JB)*FLUX(J,JB)+EMIS*ST4A(J,K,1)
620 DFABS(J,K)=DFABS(J,K)-FLUX(J,JB)
621 ENDDO
622 ENDDO
623
624 C 4.4 Outgoing longwave radiation
625
626 DO JB=1,NBAND
627 DO J=1,NGP
628 FTOP(J)=FTOP(J)+FLUX(J,JB)
629 ENDDO
630 ENDDO
631
632 DO J=1,NGP
633 FTOP(J)=FTOP(J)+OZUPP(J)
634 ENDDO
635
636 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
637 #endif /* ALLOW_AIM */
638
639 RETURN
640 END
641
642
643 SUBROUTINE RADSET( myThid )
644
645 C--
646 C-- SUBROUTINE RADSET
647 C--
648 C-- Purpose: compute energy fractions in LW bands
649 C-- as a function of temperature
650 C-- Initialized common blocks: RADFIX
651
652 IMPLICIT NONE
653
654 C Resolution parameters
655
656 C-- size for MITgcm & Physics package :
657 #include "AIM_SIZE.h"
658
659 #include "EEPARAMS.h"
660
661 C Radiation constants
662 #include "com_radcon.h"
663
664 C-- Routine arguments:
665 INTEGER myThid
666
667 #ifdef ALLOW_AIM
668
669 C-- Local variables:
670 INTEGER JTEMP, JB
671 _RL EPS3
672
673 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
674
675 EPS3=0.95 _d 0
676
677 DO JTEMP=200,320
678 FBAND(JTEMP,0)= EPSLW
679 FBAND(JTEMP,2)= 0.148 _d 0 - 3.0 _d -6 *(JTEMP-247)**2
680 FBAND(JTEMP,3)=(0.375 _d 0 - 5.5 _d -6 *(JTEMP-282)**2)*EPS3
681 FBAND(JTEMP,4)= 0.314 _d 0 + 1.0 _d -5 *(JTEMP-315)**2
682 FBAND(JTEMP,1)= 1. _d 0 -(FBAND(JTEMP,0)+FBAND(JTEMP,2)
683 & +FBAND(JTEMP,3)+FBAND(JTEMP,4))
684 ENDDO
685
686 DO JB=0,NBAND
687 DO JTEMP=lwTemp1,199
688 FBAND(JTEMP,JB)=FBAND(200,JB)
689 ENDDO
690 DO JTEMP=321,lwTemp2
691 FBAND(JTEMP,JB)=FBAND(320,JB)
692 ENDDO
693 ENDDO
694
695 C---+----1----+----2----+----3----+----4----+----5----+----6----+----7-|--+----|
696 #endif /* ALLOW_AIM */
697
698 RETURN
699 END

  ViewVC Help
Powered by ViewVC 1.1.22