Line data Source code
1 : ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_lw/src/rrtmg_lw_rtrnmc.f90,v $
2 : ! author: $Author: mike $
3 : ! revision: $Revision: 1.3 $
4 : ! created: $Date: 2008/04/24 16:17:28 $
5 : !
6 : module rrtmg_lw_rtrnmc
7 :
8 : ! --------------------------------------------------------------------------
9 : ! | |
10 : ! | Copyright 2002-2007, Atmospheric & Environmental Research, Inc. (AER). |
11 : ! | This software may be used, copied, or redistributed as long as it is |
12 : ! | not sold and this copyright notice is reproduced on each copy made. |
13 : ! | This model is provided as is without any express or implied warranties. |
14 : ! | (http://www.rtweb.aer.com/) |
15 : ! | |
16 : ! --------------------------------------------------------------------------
17 :
18 : ! --------- Modules ----------
19 :
20 : use shr_kind_mod, only: r8 => shr_kind_r8
21 :
22 : use parrrtm, only: mg, nbndlw, ngptlw
23 : use rrlw_con, only: fluxfac, heatfac
24 : use rrlw_wvn, only: delwave, ngb, ngs
25 : use rrlw_tbl, only: tblint, bpade, tau_tbl, exp_tbl, tfn_tbl
26 :
27 : implicit none
28 :
29 : contains
30 :
31 : !-----------------------------------------------------------------------------
32 55296 : subroutine rtrnmc(nlayers, istart, iend, iout, pz, semiss, &
33 27648 : cldfmc, taucmc, planklay, planklev, plankbnd, &
34 27648 : pwvcm, fracs, taut, &
35 27648 : totuflux, totdflux, fnet, htr, &
36 27648 : totuclfl, totdclfl, fnetc, htrc, totufluxs, totdfluxs )
37 : !-----------------------------------------------------------------------------
38 : !
39 : ! Original version: E. J. Mlawer, et al. RRTM_V3.0
40 : ! Revision for GCMs: Michael J. Iacono; October, 2002
41 : ! Revision for F90: Michael J. Iacono; June, 2006
42 : !
43 : ! This program calculates the upward fluxes, downward fluxes, and
44 : ! heating rates for an arbitrary clear or cloudy atmosphere. The input
45 : ! to this program is the atmospheric profile, all Planck function
46 : ! information, and the cloud fraction by layer. A variable diffusivity
47 : ! angle (SECDIFF) is used for the angle integration. Bands 2-3 and 5-9
48 : ! use a value for SECDIFF that varies from 1.50 to 1.80 as a function of
49 : ! the column water vapor, and other bands use a value of 1.66. The Gaussian
50 : ! weight appropriate to this angle (WTDIFF=0.5) is applied here. Note that
51 : ! use of the emissivity angle for the flux integration can cause errors of
52 : ! 1 to 4 W/m2 within cloudy layers.
53 : ! Clouds are treated with the McICA stochastic approach and maximum-random
54 : ! cloud overlap.
55 : !***************************************************************************
56 :
57 : ! ------- Declarations -------
58 :
59 : ! ----- Input -----
60 : integer, intent(in) :: nlayers ! total number of layers
61 : integer, intent(in) :: istart ! beginning band of calculation
62 : integer, intent(in) :: iend ! ending band of calculation
63 : integer, intent(in) :: iout ! output option flag
64 :
65 : ! Atmosphere
66 : real(kind=r8), intent(in) :: pz(0:nlayers) ! level (interface) pressures (hPa, mb)
67 : ! Dimensions: (0:nlayers)
68 : real(kind=r8), intent(in) :: pwvcm ! precipitable water vapor (cm)
69 : real(kind=r8), intent(in) :: semiss(nbndlw) ! lw surface emissivity
70 : ! Dimensions: (nbndlw)
71 : real(kind=r8), intent(in) :: planklay(nlayers,nbndlw) !
72 : ! Dimensions: (nlayers,nbndlw)
73 : real(kind=r8), intent(in) :: planklev(0:nlayers,nbndlw) !
74 : ! Dimensions: (0:nlayers,nbndlw)
75 : real(kind=r8), intent(in) :: plankbnd(nbndlw) !
76 : ! Dimensions: (nbndlw)
77 : real(kind=r8), intent(in) :: fracs(nlayers,ngptlw) !
78 : ! Dimensions: (nlayers,ngptw)
79 : real(kind=r8), intent(in) :: taut(nlayers,ngptlw) ! gaseous + aerosol optical depths
80 : ! Dimensions: (nlayers,ngptlw)
81 :
82 : ! Clouds
83 : real(kind=r8), intent(in) :: cldfmc(ngptlw,nlayers) ! layer cloud fraction [mcica]
84 : ! Dimensions: (ngptlw,nlayers)
85 : real(kind=r8), intent(in) :: taucmc(ngptlw,nlayers) ! layer cloud optical depth [mcica]
86 : ! Dimensions: (ngptlw,nlayers)
87 :
88 : ! ----- Output -----
89 : real(kind=r8), intent(out) :: totuflux(0:) ! upward longwave flux (w/m2)
90 : ! Dimensions: (0:nlayers)
91 : real(kind=r8), intent(out) :: totdflux(0:) ! downward longwave flux (w/m2)
92 : ! Dimensions: (0:nlayers)
93 : real(kind=r8), intent(out) :: fnet(0:) ! net longwave flux (w/m2)
94 : ! Dimensions: (0:nlayers)
95 : real(kind=r8), intent(out) :: htr(0:) ! longwave heating rate (k/day)
96 : ! Dimensions: (0:nlayers)
97 : real(kind=r8), intent(out) :: totuclfl(0:) ! clear sky upward longwave flux (w/m2)
98 : ! Dimensions: (0:nlayers)
99 : real(kind=r8), intent(out) :: totdclfl(0:) ! clear sky downward longwave flux (w/m2)
100 : ! Dimensions: (0:nlayers)
101 : real(kind=r8), intent(out) :: fnetc(0:) ! clear sky net longwave flux (w/m2)
102 : ! Dimensions: (0:nlayers)
103 : real(kind=r8), intent(out) :: htrc(0:) ! clear sky longwave heating rate (k/day)
104 : ! Dimensions: (0:nlayers)
105 : real(kind=r8), intent(out) :: totufluxs(:,0:) ! upward longwave flux spectral (w/m2)
106 : ! Dimensions: (nbndlw, 0:nlayers)
107 : real(kind=r8), intent(out) :: totdfluxs(:,0:) ! downward longwave flux spectral (w/m2)
108 : ! Dimensions: (nbndlw, 0:nlayers)
109 :
110 : ! ----- Local -----
111 : ! Declarations for radiative transfer
112 55296 : real(kind=r8) :: abscld(nlayers,ngptlw)
113 55296 : real(kind=r8) :: atot(nlayers)
114 55296 : real(kind=r8) :: atrans(nlayers)
115 55296 : real(kind=r8) :: bbugas(nlayers)
116 55296 : real(kind=r8) :: bbutot(nlayers)
117 55296 : real(kind=r8) :: clrurad(0:nlayers)
118 55296 : real(kind=r8) :: clrdrad(0:nlayers)
119 55296 : real(kind=r8) :: efclfrac(nlayers,ngptlw)
120 55296 : real(kind=r8) :: uflux(0:nlayers)
121 55296 : real(kind=r8) :: dflux(0:nlayers)
122 55296 : real(kind=r8) :: urad(0:nlayers)
123 55296 : real(kind=r8) :: drad(0:nlayers)
124 55296 : real(kind=r8) :: uclfl(0:nlayers)
125 55296 : real(kind=r8) :: dclfl(0:nlayers)
126 55296 : real(kind=r8) :: odcld(nlayers,ngptlw)
127 :
128 :
129 : real(kind=r8) :: secdiff(nbndlw) ! secant of diffusivity angle
130 : real(kind=r8) :: a0(nbndlw),a1(nbndlw),a2(nbndlw) ! diffusivity angle adjustment coefficients
131 : real(kind=r8) :: wtdiff, rec_6
132 : real(kind=r8) :: transcld, radld, radclrd, plfrac, blay, dplankup, dplankdn
133 : real(kind=r8) :: odepth, odtot, odepth_rec, odtot_rec, gassrc
134 : real(kind=r8) :: tblind, tfactot, bbd, bbdtot, tfacgas, transc, tausfac
135 : real(kind=r8) :: rad0, reflect, radlu, radclru
136 :
137 55296 : integer :: icldlyr(nlayers) ! flag for cloud in layer
138 : integer :: ibnd, ib, iband, lay, lev, l, ig ! loop indices
139 : integer :: igc ! g-point interval counter
140 : integer :: iclddn ! flag for cloud in down path
141 : integer :: ittot, itgas, itr ! lookup table indices
142 :
143 : ! ------- Definitions -------
144 : ! input
145 : ! nlayers ! number of model layers
146 : ! ngptlw ! total number of g-point subintervals
147 : ! nbndlw ! number of longwave spectral bands
148 : ! secdiff ! diffusivity angle
149 : ! wtdiff ! weight for radiance to flux conversion
150 : ! pavel ! layer pressures (mb)
151 : ! pz ! level (interface) pressures (mb)
152 : ! tavel ! layer temperatures (k)
153 : ! tz ! level (interface) temperatures(mb)
154 : ! tbound ! surface temperature (k)
155 : ! cldfrac ! layer cloud fraction
156 : ! taucloud ! layer cloud optical depth
157 : ! itr ! integer look-up table index
158 : ! icldlyr ! flag for cloudy layers
159 : ! iclddn ! flag for cloud in column at any layer
160 : ! semiss ! surface emissivities for each band
161 : ! reflect ! surface reflectance
162 : ! bpade ! 1/(pade constant)
163 : ! tau_tbl ! clear sky optical depth look-up table
164 : ! exp_tbl ! exponential look-up table for transmittance
165 : ! tfn_tbl ! tau transition function look-up table
166 :
167 : ! local
168 : ! atrans ! gaseous absorptivity
169 : ! abscld ! cloud absorptivity
170 : ! atot ! combined gaseous and cloud absorptivity
171 : ! odclr ! clear sky (gaseous) optical depth
172 : ! odcld ! cloud optical depth
173 : ! odtot ! optical depth of gas and cloud
174 : ! tfacgas ! gas-only pade factor, used for planck fn
175 : ! tfactot ! gas and cloud pade factor, used for planck fn
176 : ! bbdgas ! gas-only planck function for downward rt
177 : ! bbugas ! gas-only planck function for upward rt
178 : ! bbdtot ! gas and cloud planck function for downward rt
179 : ! bbutot ! gas and cloud planck function for upward calc.
180 : ! gassrc ! source radiance due to gas only
181 : ! efclfrac ! effective cloud fraction
182 : ! radlu ! spectrally summed upward radiance
183 : ! radclru ! spectrally summed clear sky upward radiance
184 : ! urad ! upward radiance by layer
185 : ! clrurad ! clear sky upward radiance by layer
186 : ! radld ! spectrally summed downward radiance
187 : ! radclrd ! spectrally summed clear sky downward radiance
188 : ! drad ! downward radiance by layer
189 : ! clrdrad ! clear sky downward radiance by layer
190 :
191 : ! output
192 : ! totuflux ! upward longwave flux (w/m2)
193 : ! totdflux ! downward longwave flux (w/m2)
194 : ! fnet ! net longwave flux (w/m2)
195 : ! htr ! longwave heating rate (k/day)
196 : ! totuclfl ! clear sky upward longwave flux (w/m2)
197 : ! totdclfl ! clear sky downward longwave flux (w/m2)
198 : ! fnetc ! clear sky net longwave flux (w/m2)
199 : ! htrc ! clear sky longwave heating rate (k/day)
200 :
201 :
202 : ! This secant and weight corresponds to the standard diffusivity
203 : ! angle. This initial value is redefined below for some bands.
204 : data wtdiff /0.5_r8/
205 : data rec_6 /0.166667_r8/
206 :
207 : ! Reset diffusivity angle for Bands 2-3 and 5-9 to vary (between 1.50
208 : ! and 1.80) as a function of total column water vapor. The function
209 : ! has been defined to minimize flux and cooling rate errors in these bands
210 : ! over a wide range of precipitable water values.
211 : data a0 / 1.66_r8, 1.55_r8, 1.58_r8, 1.66_r8, &
212 : 1.54_r8, 1.454_r8, 1.89_r8, 1.33_r8, &
213 : 1.668_r8, 1.66_r8, 1.66_r8, 1.66_r8, &
214 : 1.66_r8, 1.66_r8, 1.66_r8, 1.66_r8 /
215 : data a1 / 0.00_r8, 0.25_r8, 0.22_r8, 0.00_r8, &
216 : 0.13_r8, 0.446_r8, -0.10_r8, 0.40_r8, &
217 : -0.006_r8, 0.00_r8, 0.00_r8, 0.00_r8, &
218 : 0.00_r8, 0.00_r8, 0.00_r8, 0.00_r8 /
219 : data a2 / 0.00_r8, -12.0_r8, -11.7_r8, 0.00_r8, &
220 : -0.72_r8,-0.243_r8, 0.19_r8,-0.062_r8, &
221 : 0.414_r8, 0.00_r8, 0.00_r8, 0.00_r8, &
222 : 0.00_r8, 0.00_r8, 0.00_r8, 0.00_r8 /
223 :
224 470016 : do ibnd = 1,nbndlw
225 470016 : if (ibnd.eq.1 .or. ibnd.eq.4 .or. ibnd.ge.10) then
226 248832 : secdiff(ibnd) = 1.66_r8
227 : else
228 193536 : secdiff(ibnd) = a0(ibnd) + a1(ibnd)*exp(a2(ibnd)*pwvcm)
229 193536 : if (secdiff(ibnd) .gt. 1.80_r8) secdiff(ibnd) = 1.80_r8
230 193536 : if (secdiff(ibnd) .lt. 1.50_r8) secdiff(ibnd) = 1.50_r8
231 : endif
232 : enddo
233 :
234 27648 : urad(0) = 0.0_r8
235 27648 : drad(0) = 0.0_r8
236 27648 : totuflux(0) = 0.0_r8
237 27648 : totdflux(0) = 0.0_r8
238 27648 : clrurad(0) = 0.0_r8
239 27648 : clrdrad(0) = 0.0_r8
240 27648 : totuclfl(0) = 0.0_r8
241 27648 : totdclfl(0) = 0.0_r8
242 :
243 1741824 : do lay = 1, nlayers
244 1714176 : urad(lay) = 0.0_r8
245 1714176 : drad(lay) = 0.0_r8
246 1714176 : totuflux(lay) = 0.0_r8
247 1714176 : totdflux(lay) = 0.0_r8
248 1714176 : clrurad(lay) = 0.0_r8
249 1714176 : clrdrad(lay) = 0.0_r8
250 1714176 : totuclfl(lay) = 0.0_r8
251 1714176 : totdclfl(lay) = 0.0_r8
252 1741824 : icldlyr(lay) = 0
253 : enddo
254 : ! Change to band loop?
255 3898368 : do ig = 1, ngptlw
256 243883008 : do lay = 1, nlayers
257 243855360 : if (cldfmc(ig,lay) .eq. 1._r8) then
258 17969071 : ib = ngb(ig)
259 17969071 : odcld(lay,ig) = secdiff(ib) * taucmc(ig,lay)
260 17969071 : transcld = exp(-odcld(lay,ig))
261 17969071 : abscld(lay,ig) = 1._r8 - transcld
262 17969071 : efclfrac(lay,ig) = abscld(lay,ig) * cldfmc(ig,lay)
263 17969071 : icldlyr(lay) = 1
264 : else
265 222015569 : odcld(lay,ig) = 0.0_r8
266 222015569 : abscld(lay,ig) = 0.0_r8
267 222015569 : efclfrac(lay,ig) = 0.0_r8
268 : endif
269 : enddo
270 :
271 : enddo
272 :
273 27648 : igc = 1
274 : ! Loop over frequency bands.
275 470016 : do iband = istart, iend
276 :
277 : ! Reinitialize g-point counter for each band if output for each band is requested.
278 442368 : if (iout.gt.0.and.iband.ge.2) igc = ngs(iband-1)+1
279 :
280 : ! Loop over g-channels.
281 : 1000 continue
282 :
283 : ! Radiative transfer starts here.
284 3870720 : radld = 0._r8
285 3870720 : radclrd = 0._r8
286 3870720 : iclddn = 0
287 :
288 : ! Downward radiative transfer loop.
289 :
290 243855360 : do lev = nlayers, 1, -1
291 239984640 : plfrac = fracs(lev,igc)
292 239984640 : blay = planklay(lev,iband)
293 239984640 : dplankup = planklev(lev,iband) - blay
294 239984640 : dplankdn = planklev(lev-1,iband) - blay
295 239984640 : odepth = secdiff(iband) * taut(lev,igc)
296 239984640 : if (odepth .lt. 0.0_r8) odepth = 0.0_r8
297 : ! Cloudy layer
298 239984640 : if (icldlyr(lev).eq.1) then
299 33346600 : iclddn = 1
300 33346600 : odtot = odepth + odcld(lev,igc)
301 33346600 : if (odtot .lt. 0.06_r8) then
302 4864239 : atrans(lev) = odepth - 0.5_r8*odepth*odepth
303 4864239 : odepth_rec = rec_6*odepth
304 4864239 : gassrc = plfrac*(blay+dplankdn*odepth_rec)*atrans(lev)
305 :
306 4864239 : atot(lev) = odtot - 0.5_r8*odtot*odtot
307 4864239 : odtot_rec = rec_6*odtot
308 4864239 : bbdtot = plfrac * (blay+dplankdn*odtot_rec)
309 4864239 : bbd = plfrac*(blay+dplankdn*odepth_rec)
310 : radld = radld - radld * (atrans(lev) + &
311 : efclfrac(lev,igc) * (1. - atrans(lev))) + &
312 : gassrc + cldfmc(igc,lev) * &
313 4864239 : (bbdtot * atot(lev) - gassrc)
314 4864239 : drad(lev-1) = drad(lev-1) + radld
315 :
316 4864239 : bbugas(lev) = plfrac * (blay+dplankup*odepth_rec)
317 4864239 : bbutot(lev) = plfrac * (blay+dplankup*odtot_rec)
318 :
319 28482361 : elseif (odepth .le. 0.06_r8) then
320 3087893 : atrans(lev) = odepth - 0.5_r8*odepth*odepth
321 3087893 : odepth_rec = rec_6*odepth
322 3087893 : gassrc = plfrac*(blay+dplankdn*odepth_rec)*atrans(lev)
323 :
324 3087893 : odtot = odepth + odcld(lev,igc)
325 3087893 : tblind = odtot/(bpade+odtot)
326 3087893 : ittot = tblint*tblind + 0.5_r8
327 3087893 : tfactot = tfn_tbl(ittot)
328 3087893 : bbdtot = plfrac * (blay + tfactot*dplankdn)
329 3087893 : bbd = plfrac*(blay+dplankdn*odepth_rec)
330 3087893 : atot(lev) = 1. - exp_tbl(ittot)
331 :
332 : radld = radld - radld * (atrans(lev) + &
333 : efclfrac(lev,igc) * (1._r8 - atrans(lev))) + &
334 : gassrc + cldfmc(igc,lev) * &
335 3087893 : (bbdtot * atot(lev) - gassrc)
336 3087893 : drad(lev-1) = drad(lev-1) + radld
337 :
338 3087893 : bbugas(lev) = plfrac * (blay + dplankup*odepth_rec)
339 3087893 : bbutot(lev) = plfrac * (blay + tfactot * dplankup)
340 :
341 : else
342 :
343 25394468 : tblind = odepth/(bpade+odepth)
344 25394468 : itgas = tblint*tblind+0.5_r8
345 25394468 : odepth = tau_tbl(itgas)
346 25394468 : atrans(lev) = 1._r8 - exp_tbl(itgas)
347 25394468 : tfacgas = tfn_tbl(itgas)
348 25394468 : gassrc = atrans(lev) * plfrac * (blay + tfacgas*dplankdn)
349 :
350 25394468 : odtot = odepth + odcld(lev,igc)
351 25394468 : tblind = odtot/(bpade+odtot)
352 25394468 : ittot = tblint*tblind + 0.5_r8
353 25394468 : tfactot = tfn_tbl(ittot)
354 25394468 : bbdtot = plfrac * (blay + tfactot*dplankdn)
355 25394468 : bbd = plfrac*(blay+tfacgas*dplankdn)
356 25394468 : atot(lev) = 1._r8 - exp_tbl(ittot)
357 :
358 : radld = radld - radld * (atrans(lev) + &
359 : efclfrac(lev,igc) * (1._r8 - atrans(lev))) + &
360 : gassrc + cldfmc(igc,lev) * &
361 25394468 : (bbdtot * atot(lev) - gassrc)
362 25394468 : drad(lev-1) = drad(lev-1) + radld
363 25394468 : bbugas(lev) = plfrac * (blay + tfacgas * dplankup)
364 25394468 : bbutot(lev) = plfrac * (blay + tfactot * dplankup)
365 : endif
366 : ! Clear layer
367 : else
368 206638040 : if (odepth .le. 0.06_r8) then
369 140085174 : atrans(lev) = odepth-0.5_r8*odepth*odepth
370 140085174 : odepth = rec_6*odepth
371 140085174 : bbd = plfrac*(blay+dplankdn*odepth)
372 140085174 : bbugas(lev) = plfrac*(blay+dplankup*odepth)
373 : else
374 66552866 : tblind = odepth/(bpade+odepth)
375 66552866 : itr = tblint*tblind+0.5_r8
376 66552866 : transc = exp_tbl(itr)
377 66552866 : atrans(lev) = 1._r8-transc
378 66552866 : tausfac = tfn_tbl(itr)
379 66552866 : bbd = plfrac*(blay+tausfac*dplankdn)
380 66552866 : bbugas(lev) = plfrac * (blay + tausfac * dplankup)
381 : endif
382 206638040 : radld = radld + (bbd-radld)*atrans(lev)
383 206638040 : drad(lev-1) = drad(lev-1) + radld
384 : endif
385 : ! Set clear sky stream to total sky stream as long as layers
386 : ! remain clear. Streams diverge when a cloud is reached (iclddn=1),
387 : ! and clear sky stream must be computed separately from that point.
388 243855360 : if (iclddn.eq.1) then
389 48095600 : radclrd = radclrd + (bbd-radclrd) * atrans(lev)
390 48095600 : clrdrad(lev-1) = clrdrad(lev-1) + radclrd
391 : else
392 191889040 : radclrd = radld
393 191889040 : clrdrad(lev-1) = drad(lev-1)
394 : endif
395 : enddo
396 :
397 : ! Spectral emissivity & reflectance
398 : ! Include the contribution of spectrally varying longwave emissivity
399 : ! and reflection from the surface to the upward radiative transfer.
400 : ! Note: Spectral and Lambertian reflection are identical for the
401 : ! diffusivity angle flux integration used here.
402 :
403 3870720 : rad0 = fracs(1,igc) * plankbnd(iband)
404 : ! Add in specular reflection of surface downward radiance.
405 3870720 : reflect = 1._r8 - semiss(iband)
406 3870720 : radlu = rad0 + reflect * radld
407 3870720 : radclru = rad0 + reflect * radclrd
408 :
409 :
410 : ! Upward radiative transfer loop.
411 3870720 : urad(0) = urad(0) + radlu
412 3870720 : clrurad(0) = clrurad(0) + radclru
413 :
414 243855360 : do lev = 1, nlayers
415 : ! Cloudy layer
416 239984640 : if (icldlyr(lev) .eq. 1) then
417 33346600 : gassrc = bbugas(lev) * atrans(lev)
418 : radlu = radlu - radlu * (atrans(lev) + &
419 : efclfrac(lev,igc) * (1._r8 - atrans(lev))) + &
420 : gassrc + cldfmc(igc,lev) * &
421 33346600 : (bbutot(lev) * atot(lev) - gassrc)
422 33346600 : urad(lev) = urad(lev) + radlu
423 : ! Clear layer
424 : else
425 206638040 : radlu = radlu + (bbugas(lev)-radlu)*atrans(lev)
426 206638040 : urad(lev) = urad(lev) + radlu
427 : endif
428 : ! Set clear sky stream to total sky stream as long as all layers
429 : ! are clear (iclddn=0). Streams must be calculated separately at
430 : ! all layers when a cloud is present (ICLDDN=1), because surface
431 : ! reflectance is different for each stream.
432 243855360 : if (iclddn.eq.1) then
433 216800360 : radclru = radclru + (bbugas(lev)-radclru)*atrans(lev)
434 216800360 : clrurad(lev) = clrurad(lev) + radclru
435 : else
436 23184280 : radclru = radlu
437 23184280 : clrurad(lev) = urad(lev)
438 : endif
439 : enddo
440 :
441 : ! Increment g-point counter
442 3870720 : igc = igc + 1
443 : ! Return to continue radiative transfer for all g-channels in present band
444 3870720 : if (igc .le. ngs(iband)) go to 1000
445 :
446 : ! Process longwave output from band for total and clear streams.
447 : ! Calculate upward, downward, and net flux.
448 28311552 : do lev = nlayers, 0, -1
449 27869184 : uflux(lev) = urad(lev)*wtdiff
450 27869184 : dflux(lev) = drad(lev)*wtdiff
451 27869184 : urad(lev) = 0.0_r8
452 27869184 : drad(lev) = 0.0_r8
453 27869184 : uclfl(lev) = clrurad(lev)*wtdiff
454 27869184 : dclfl(lev) = clrdrad(lev)*wtdiff
455 27869184 : clrurad(lev) = 0.0_r8
456 28311552 : clrdrad(lev) = 0.0_r8
457 : enddo
458 :
459 28339200 : do lev = nlayers, 0, -1
460 27869184 : totuflux(lev) = totuflux(lev) + uflux(lev) * delwave(iband)
461 27869184 : totdflux(lev) = totdflux(lev) + dflux(lev) * delwave(iband)
462 27869184 : totuclfl(lev) = totuclfl(lev) + uclfl(lev) * delwave(iband)
463 27869184 : totdclfl(lev) = totdclfl(lev) + dclfl(lev) * delwave(iband)
464 27869184 : totufluxs(iband,lev) = uflux(lev) * delwave(iband)
465 28311552 : totdfluxs(iband,lev) = dflux(lev) * delwave(iband)
466 : enddo
467 : ! End spectral band loop
468 : enddo
469 :
470 : ! Calculate fluxes at surface
471 27648 : totuflux(0) = totuflux(0) * fluxfac
472 27648 : totdflux(0) = totdflux(0) * fluxfac
473 470016 : totufluxs(:,0) = totufluxs(:,0) * fluxfac
474 470016 : totdfluxs(:,0) = totdfluxs(:,0) * fluxfac
475 27648 : fnet(0) = totuflux(0) - totdflux(0)
476 27648 : totuclfl(0) = totuclfl(0) * fluxfac
477 27648 : totdclfl(0) = totdclfl(0) * fluxfac
478 27648 : fnetc(0) = totuclfl(0) - totdclfl(0)
479 :
480 : ! Calculate fluxes at model levels
481 1741824 : do lev = 1, nlayers
482 1714176 : totuflux(lev) = totuflux(lev) * fluxfac
483 1714176 : totdflux(lev) = totdflux(lev) * fluxfac
484 29140992 : totufluxs(:,lev) = totufluxs(:,lev) * fluxfac
485 29140992 : totdfluxs(:,lev) = totdfluxs(:,lev) * fluxfac
486 1714176 : fnet(lev) = totuflux(lev) - totdflux(lev)
487 1714176 : totuclfl(lev) = totuclfl(lev) * fluxfac
488 1714176 : totdclfl(lev) = totdclfl(lev) * fluxfac
489 1714176 : fnetc(lev) = totuclfl(lev) - totdclfl(lev)
490 1714176 : l = lev - 1
491 :
492 : ! Calculate heating rates at model layers
493 1714176 : htr(l)=heatfac*(fnet(l)-fnet(lev))/(pz(l)-pz(lev))
494 1741824 : htrc(l)=heatfac*(fnetc(l)-fnetc(lev))/(pz(l)-pz(lev))
495 : enddo
496 :
497 : ! Set heating rate to zero in top layer
498 27648 : htr(nlayers) = 0.0_r8
499 27648 : htrc(nlayers) = 0.0_r8
500 :
501 27648 : end subroutine rtrnmc
502 :
503 : end module rrtmg_lw_rtrnmc
504 :
|