Line data Source code
1 : ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_lw/src/rrtmg_lw.f90,v $
2 : ! author: $Author: mike $
3 : ! revision: $Revision: 1.6 $
4 : ! created: $Date: 2008/04/24 16:17:27 $
5 : !
6 :
7 : module rrtmg_lw_rad
8 :
9 : ! --------------------------------------------------------------------------
10 : ! | |
11 : ! | Copyright 2002-2007, Atmospheric & Environmental Research, Inc. (AER). |
12 : ! | This software may be used, copied, or redistributed as long as it is |
13 : ! | not sold and this copyright notice is reproduced on each copy made. |
14 : ! | This model is provided as is without any express or implied warranties. |
15 : ! | (http://www.rtweb.aer.com/) |
16 : ! | |
17 : ! --------------------------------------------------------------------------
18 : !
19 : ! ****************************************************************************
20 : ! * *
21 : ! * RRTMG_LW *
22 : ! * *
23 : ! * *
24 : ! * *
25 : ! * a rapid radiative transfer model *
26 : ! * for the longwave region *
27 : ! * for application to general circulation models *
28 : ! * *
29 : ! * *
30 : ! * Atmospheric and Environmental Research, Inc. *
31 : ! * 131 Hartwell Avenue *
32 : ! * Lexington, MA 02421 *
33 : ! * *
34 : ! * *
35 : ! * Eli J. Mlawer *
36 : ! * Jennifer S. Delamere *
37 : ! * Michael J. Iacono *
38 : ! * Shepard A. Clough *
39 : ! * *
40 : ! * *
41 : ! * *
42 : ! * *
43 : ! * *
44 : ! * *
45 : ! * email: miacono@aer.com *
46 : ! * email: emlawer@aer.com *
47 : ! * email: jdelamer@aer.com *
48 : ! * *
49 : ! * The authors wish to acknowledge the contributions of the *
50 : ! * following people: Steven J. Taubman, Karen Cady-Pereira, *
51 : ! * Patrick D. Brown, Ronald E. Farren, Luke Chen, Robert Bergstrom. *
52 : ! * *
53 : ! ****************************************************************************
54 :
55 : ! -------- Modules --------
56 :
57 : use shr_kind_mod, only: r8=>shr_kind_r8
58 :
59 : use mcica_subcol_gen_lw, only: mcica_subcol_lw
60 : use rrtmg_lw_setcoef, only: setcoef
61 : use rrtmg_lw_taumol, only: taumol
62 : use rrtmg_lw_rtrnmc, only: rtrnmc
63 :
64 : implicit none
65 :
66 : public :: rrtmg_lw
67 :
68 : ! Set iaer to select aerosol option
69 : ! iaer = 0, no aerosols
70 : ! iaer = 10, input total aerosol optical depth (tauaer) directly
71 : integer, parameter :: iaer = 10
72 :
73 : !=========================================================================================
74 : contains
75 : !=========================================================================================
76 :
77 2304 : subroutine rrtmg_lw &
78 : (lchnk ,ncol ,nlay ,icld , &
79 2304 : play ,plev ,tlay ,tlev ,tsfc ,h2ovmr , &
80 2304 : o3vmr ,co2vmr ,ch4vmr ,o2vmr ,n2ovmr ,&
81 2304 : cfc11vmr,cfc12vmr, &
82 2304 : cfc22vmr,ccl4vmr ,emis , &
83 2304 : cldfmcl ,taucmcl ,ciwpmcl ,clwpmcl ,reicmcl ,relqmcl , &
84 2304 : tauaer , &
85 2304 : uflx ,dflx ,hr ,uflxc ,dflxc, hrc, uflxs, dflxs )
86 :
87 : ! -------- Description --------
88 :
89 : ! This program is the driver subroutine for RRTMG_LW, the AER LW radiation
90 : ! model for application to GCMs, that has been adapted from RRTM_LW for
91 : ! improved efficiency.
92 : !
93 : ! This routine:
94 : ! a) calls INATM to read in the atmospheric profile from GCM;
95 : ! all layering in RRTMG is ordered from surface to toa.
96 : ! b) call to CLDPRMC removed -- CAM supplies cloud optical depths
97 : ! c) calls SETCOEF to calculate various quantities needed for
98 : ! the radiative transfer algorithm
99 : ! d) calls TAUMOL to calculate gaseous optical depths for each
100 : ! of the 16 spectral bands
101 : ! e) calls RTRNMC (for both clear and cloudy profiles) to perform the
102 : ! radiative transfer calculation using McICA, the Monte-Carlo
103 : ! Independent Column Approximation, to represent sub-grid scale
104 : ! cloud variability
105 : ! f) passes the necessary fluxes and cooling rates back to GCM
106 : !
107 : ! *** This version uses McICA ***
108 : ! Monte Carlo Independent Column Approximation (McICA, Pincus et al.,
109 : ! JC, 2003) method is applied to the forward model calculation
110 : !
111 : ! This call to RRTMG_LW must be preceeded by a call to the module
112 : ! mcica_subcol_gen_lw.f90 to run the McICA sub-column cloud generator,
113 : ! which will provide the cloud physical or cloud optical properties
114 : ! on the RRTMG quadrature point (ngpt) dimension.
115 : !
116 : ! *** This version requires that cloud optical properties be input ***
117 : !
118 : ! *** This version requires that aerosol optical properties be input ***
119 : ! Input aerosol optical depth directly by layer and spectral band (iaer=10);
120 : ! band average optical depth at the mid-point of each spectral band.
121 : ! RRTMG_LW currently treats only aerosol absorption;
122 : ! scattering capability is not presently available.
123 : !
124 : !
125 : ! ------- Modifications -------
126 : !
127 : ! This version of RRTMG_LW has been modified from RRTM_LW to use a reduced
128 : ! set of g-points for application to GCMs.
129 : !
130 : !-- Original version (derived from RRTM_LW), reduction of g-points, other
131 : ! revisions for use with GCMs.
132 : ! 1999: M. J. Iacono, AER, Inc.
133 : !-- Adapted for use with NCAR/CAM.
134 : ! May 2004: M. J. Iacono, AER, Inc.
135 : !-- Revised to add McICA capability.
136 : ! Nov 2005: M. J. Iacono, AER, Inc.
137 : !-- Conversion to F90 formatting for consistency with rrtmg_sw.
138 : ! Feb 2007: M. J. Iacono, AER, Inc.
139 : !-- Modifications to formatting to use assumed-shape arrays.
140 : ! Aug 2007: M. J. Iacono, AER, Inc.
141 : !-- Modified to add longwave aerosol absorption.
142 : ! Apr 2008: M. J. Iacono, AER, Inc.
143 :
144 : use parrrtm, only: nbndlw, ngptlw, maxxsec, mxmol
145 : use rrlw_con, only: fluxfac, oneminus, pi
146 : use rrlw_wvn, only: ngb
147 :
148 : ! ----- Input -----
149 : integer, intent(in) :: lchnk ! chunk identifier
150 : integer, intent(in) :: ncol ! Number of horizontal columns
151 : integer, intent(in) :: nlay ! Number of model layers
152 : integer, intent(in) :: icld ! Cloud overlap method
153 : ! 0: Clear only
154 : ! 1: Random
155 : ! 2: Maximum/random
156 : ! 3: Maximum
157 : real(kind=r8), intent(in) :: play(:,:) ! Layer pressures (hPa, mb)
158 : ! Dimensions: (ncol,nlay)
159 : real(kind=r8), intent(in) :: plev(:,:) ! Interface pressures (hPa, mb)
160 : ! Dimensions: (ncol,nlay+1)
161 : real(kind=r8), intent(in) :: tlay(:,:) ! Layer temperatures (K)
162 : ! Dimensions: (ncol,nlay)
163 : real(kind=r8), intent(in) :: tlev(:,:) ! Interface temperatures (K)
164 : ! Dimensions: (ncol,nlay+1)
165 : real(kind=r8), intent(in) :: tsfc(:) ! Surface temperature (K)
166 : ! Dimensions: (ncol)
167 : real(kind=r8), intent(in) :: h2ovmr(:,:) ! H2O volume mixing ratio
168 : ! Dimensions: (ncol,nlay)
169 : real(kind=r8), intent(in) :: o3vmr(:,:) ! O3 volume mixing ratio
170 : ! Dimensions: (ncol,nlay)
171 : real(kind=r8), intent(in) :: co2vmr(:,:) ! CO2 volume mixing ratio
172 : ! Dimensions: (ncol,nlay)
173 : real(kind=r8), intent(in) :: ch4vmr(:,:) ! Methane volume mixing ratio
174 : ! Dimensions: (ncol,nlay)
175 : real(kind=r8), intent(in) :: o2vmr(:,:) ! O2 volume mixing ratio
176 : ! Dimensions: (ncol,nlay)
177 : real(kind=r8), intent(in) :: n2ovmr(:,:) ! Nitrous oxide volume mixing ratio
178 : ! Dimensions: (ncol,nlay)
179 : real(kind=r8), intent(in) :: cfc11vmr(:,:) ! CFC11 volume mixing ratio
180 : ! Dimensions: (ncol,nlay)
181 : real(kind=r8), intent(in) :: cfc12vmr(:,:) ! CFC12 volume mixing ratio
182 : ! Dimensions: (ncol,nlay)
183 : real(kind=r8), intent(in) :: cfc22vmr(:,:) ! CFC22 volume mixing ratio
184 : ! Dimensions: (ncol,nlay)
185 : real(kind=r8), intent(in) :: ccl4vmr(:,:) ! CCL4 volume mixing ratio
186 : ! Dimensions: (ncol,nlay)
187 : real(kind=r8), intent(in) :: emis(:,:) ! Surface emissivity
188 : ! Dimensions: (ncol,nbndlw)
189 :
190 : real(kind=r8), intent(in) :: cldfmcl(:,:,:) ! Cloud fraction
191 : ! Dimensions: (ngptlw,ncol,nlay)
192 : real(kind=r8), intent(in) :: ciwpmcl(:,:,:) ! Cloud ice water path (g/m2)
193 : ! Dimensions: (ngptlw,ncol,nlay)
194 : real(kind=r8), intent(in) :: clwpmcl(:,:,:) ! Cloud liquid water path (g/m2)
195 : ! Dimensions: (ngptlw,ncol,nlay)
196 : real(kind=r8), intent(in) :: reicmcl(:,:) ! Cloud ice effective radius (microns)
197 : ! Dimensions: (ncol,nlay)
198 : real(kind=r8), intent(in) :: relqmcl(:,:) ! Cloud water drop effective radius (microns)
199 : ! Dimensions: (ncol,nlay)
200 : real(kind=r8), intent(in) :: taucmcl(:,:,:) ! Cloud optical depth
201 : ! Dimensions: (ngptlw,ncol,nlay)
202 : real(kind=r8), intent(in) :: tauaer(:,:,:) ! aerosol optical depth
203 : ! at mid-point of LW spectral bands
204 : ! Dimensions: (ncol,nlay,nbndlw)
205 :
206 : ! ----- Output -----
207 :
208 : real(kind=r8), intent(out) :: uflx(:,:) ! Total sky longwave upward flux (W/m2)
209 : ! Dimensions: (ncol,nlay+1)
210 : real(kind=r8), intent(out) :: dflx(:,:) ! Total sky longwave downward flux (W/m2)
211 : ! Dimensions: (ncol,nlay+1)
212 : real(kind=r8), intent(out) :: hr(:,:) ! Total sky longwave radiative heating rate (K/d)
213 : ! Dimensions: (ncol,nlay)
214 : real(kind=r8), intent(out) :: uflxc(:,:) ! Clear sky longwave upward flux (W/m2)
215 : ! Dimensions: (ncol,nlay+1)
216 : real(kind=r8), intent(out) :: dflxc(:,:) ! Clear sky longwave downward flux (W/m2)
217 : ! Dimensions: (ncol,nlay+1)
218 : real(kind=r8), intent(out) :: hrc(:,:) ! Clear sky longwave radiative heating rate (K/d)
219 : ! Dimensions: (ncol,nlay)
220 : real(kind=r8), intent(out) :: uflxs(:,:,:) ! Total sky longwave upward flux spectral (W/m2)
221 : ! Dimensions: (nbndlw,ncol,nlay+1)
222 : real(kind=r8), intent(out) :: dflxs(:,:,:) ! Total sky longwave downward flux spectral (W/m2)
223 : ! Dimensions: (nbndlw,ncol,nlay+1)
224 :
225 : ! ----- Local -----
226 :
227 : ! Control
228 : integer :: istart ! beginning band of calculation
229 : integer :: iend ! ending band of calculation
230 : integer :: iout ! output option flag (inactive)
231 : integer :: iplon ! column loop index
232 : integer :: ims ! value for changing mcica permute seed
233 : integer :: k ! layer loop index
234 : integer :: ig ! g-point loop index
235 :
236 : ! Atmosphere
237 4608 : real(kind=r8) :: pavel(nlay) ! layer pressures (mb)
238 4608 : real(kind=r8) :: tavel(nlay) ! layer temperatures (K)
239 4608 : real(kind=r8) :: pz(0:nlay) ! level (interface) pressures (hPa, mb)
240 4608 : real(kind=r8) :: tz(0:nlay) ! level (interface) temperatures (K)
241 : real(kind=r8) :: tbound ! surface temperature (K)
242 4608 : real(kind=r8) :: coldry(nlay) ! dry air column density (mol/cm2)
243 4608 : real(kind=r8) :: wbrodl(nlay) ! broadening gas column density (mol/cm2)
244 4608 : real(kind=r8) :: wkl(mxmol,nlay) ! molecular amounts (mol/cm-2)
245 4608 : real(kind=r8) :: wx(maxxsec,nlay) ! cross-section amounts (mol/cm-2)
246 : real(kind=r8) :: pwvcm ! precipitable water vapor (cm)
247 : real(kind=r8) :: semiss(nbndlw) ! lw surface emissivity
248 4608 : real(kind=r8) :: fracs(nlay,ngptlw) !
249 4608 : real(kind=r8) :: taug(nlay,ngptlw) ! gaseous optical depths
250 4608 : real(kind=r8) :: taut(nlay,ngptlw) ! gaseous + aerosol optical depths
251 :
252 4608 : real(kind=r8) :: taua(nlay,nbndlw) ! aerosol optical depth
253 :
254 : ! Atmosphere - setcoef
255 : integer :: laytrop ! tropopause layer index
256 4608 : integer :: jp(nlay) ! lookup table index
257 4608 : integer :: jt(nlay) ! lookup table index
258 4608 : integer :: jt1(nlay) ! lookup table index
259 4608 : real(kind=r8) :: planklay(nlay,nbndlw) !
260 4608 : real(kind=r8) :: planklev(0:nlay,nbndlw) !
261 : real(kind=r8) :: plankbnd(nbndlw) !
262 :
263 4608 : real(kind=r8) :: colh2o(nlay) ! column amount (h2o)
264 4608 : real(kind=r8) :: colco2(nlay) ! column amount (co2)
265 4608 : real(kind=r8) :: colo3(nlay) ! column amount (o3)
266 4608 : real(kind=r8) :: coln2o(nlay) ! column amount (n2o)
267 4608 : real(kind=r8) :: colco(nlay) ! column amount (co)
268 4608 : real(kind=r8) :: colch4(nlay) ! column amount (ch4)
269 4608 : real(kind=r8) :: colo2(nlay) ! column amount (o2)
270 4608 : real(kind=r8) :: colbrd(nlay) ! column amount (broadening gases)
271 :
272 4608 : integer :: indself(nlay)
273 4608 : integer :: indfor(nlay)
274 4608 : real(kind=r8) :: selffac(nlay)
275 4608 : real(kind=r8) :: selffrac(nlay)
276 4608 : real(kind=r8) :: forfac(nlay)
277 4608 : real(kind=r8) :: forfrac(nlay)
278 :
279 4608 : integer :: indminor(nlay)
280 4608 : real(kind=r8) :: minorfrac(nlay)
281 4608 : real(kind=r8) :: scaleminor(nlay)
282 4608 : real(kind=r8) :: scaleminorn2(nlay)
283 :
284 : real(kind=r8) :: & !
285 4608 : fac00(nlay), fac01(nlay), &
286 4608 : fac10(nlay), fac11(nlay)
287 : real(kind=r8) :: & !
288 4608 : rat_h2oco2(nlay),rat_h2oco2_1(nlay), &
289 4608 : rat_h2oo3(nlay),rat_h2oo3_1(nlay), &
290 4608 : rat_h2on2o(nlay),rat_h2on2o_1(nlay), &
291 4608 : rat_h2och4(nlay),rat_h2och4_1(nlay), &
292 4608 : rat_n2oco2(nlay),rat_n2oco2_1(nlay), &
293 4608 : rat_o3co2(nlay),rat_o3co2_1(nlay)
294 :
295 : ! Atmosphere/clouds - cldprmc [mcica]
296 4608 : real(kind=r8) :: cldfmc(ngptlw,nlay) ! cloud fraction [mcica]
297 4608 : real(kind=r8) :: ciwpmc(ngptlw,nlay) ! cloud ice water path [mcica]
298 4608 : real(kind=r8) :: clwpmc(ngptlw,nlay) ! cloud liquid water path [mcica]
299 4608 : real(kind=r8) :: relqmc(nlay) ! liquid particle size (microns)
300 4608 : real(kind=r8) :: reicmc(nlay) ! ice particle effective radius (microns)
301 4608 : real(kind=r8) :: dgesmc(nlay) ! ice particle generalized effective size (microns)
302 4608 : real(kind=r8) :: taucmc(ngptlw,nlay) ! cloud optical depth [mcica]
303 :
304 : ! Output
305 4608 : real(kind=r8) :: totuflux(0:nlay) ! upward longwave flux (w/m2)
306 4608 : real(kind=r8) :: totdflux(0:nlay) ! downward longwave flux (w/m2)
307 4608 : real(kind=r8) :: totufluxs(nbndlw,0:nlay) ! upward longwave flux spectral (w/m2)
308 4608 : real(kind=r8) :: totdfluxs(nbndlw,0:nlay) ! downward longwave flux spectral (w/m2)
309 4608 : real(kind=r8) :: fnet(0:nlay) ! net longwave flux (w/m2)
310 4608 : real(kind=r8) :: htr(0:nlay) ! longwave heating rate (k/day)
311 4608 : real(kind=r8) :: totuclfl(0:nlay) ! clear sky upward longwave flux (w/m2)
312 4608 : real(kind=r8) :: totdclfl(0:nlay) ! clear sky downward longwave flux (w/m2)
313 4608 : real(kind=r8) :: fnetc(0:nlay) ! clear sky net longwave flux (w/m2)
314 4608 : real(kind=r8) :: htrc(0:nlay) ! clear sky longwave heating rate (k/day)
315 : !----------------------------------------------------------------------------
316 :
317 2304 : oneminus = 1._r8 - 1.e-6_r8
318 2304 : pi = 2._r8 * asin(1._r8)
319 2304 : fluxfac = pi * 2.e4_r8
320 2304 : istart = 1
321 2304 : iend = 16
322 2304 : iout = 0
323 2304 : ims = 1
324 :
325 : ! Main longitude/column loop within RRTMG.
326 29952 : do iplon = 1, ncol
327 :
328 : ! Prepare atmospheric profile from GCM for use in RRTMG, and define
329 : ! other input parameters.
330 :
331 : call inatm(iplon, nlay, icld, iaer, &
332 : play, plev, tlay, tlev, tsfc, h2ovmr, &
333 : o3vmr, co2vmr, ch4vmr, o2vmr, n2ovmr, cfc11vmr, cfc12vmr, &
334 : cfc22vmr, ccl4vmr, emis, &
335 : cldfmcl, taucmcl, ciwpmcl, clwpmcl, reicmcl, relqmcl, tauaer, &
336 : pavel, pz, tavel, tz, tbound, semiss, coldry, &
337 : wkl, wbrodl, wx, pwvcm, &
338 27648 : cldfmc, taucmc, ciwpmc, clwpmc, reicmc, dgesmc, relqmc, taua)
339 :
340 : ! Calculate information needed by the radiative transfer routine
341 : ! that is specific to this atmosphere
342 : ! by interpolating data from stored reference atmospheres.
343 :
344 : call setcoef(nlay, istart, pavel, tavel, tz, tbound, semiss, &
345 : coldry, wkl, wbrodl, &
346 : laytrop, jp, jt, jt1, planklay, planklev, plankbnd, &
347 : colh2o, colco2, colo3, coln2o, colco, colch4, colo2, &
348 : colbrd, fac00, fac01, fac10, fac11, &
349 : rat_h2oco2, rat_h2oco2_1, rat_h2oo3, rat_h2oo3_1, &
350 : rat_h2on2o, rat_h2on2o_1, rat_h2och4, rat_h2och4_1, &
351 : rat_n2oco2, rat_n2oco2_1, rat_o3co2, rat_o3co2_1, &
352 : selffac, selffrac, indself, forfac, forfrac, indfor, &
353 27648 : minorfrac, scaleminor, scaleminorn2, indminor)
354 :
355 : ! Calculate the gaseous optical depths and Planck fractions for
356 : ! each longwave spectral band.
357 :
358 : call taumol(nlay, pavel, wx, coldry, &
359 : laytrop, jp, jt, jt1, planklay, planklev, plankbnd, &
360 : colh2o, colco2, colo3, coln2o, colco, colch4, colo2, &
361 : colbrd, fac00, fac01, fac10, fac11, &
362 : rat_h2oco2, rat_h2oco2_1, rat_h2oo3, rat_h2oo3_1, &
363 : rat_h2on2o, rat_h2on2o_1, rat_h2och4, rat_h2och4_1, &
364 : rat_n2oco2, rat_n2oco2_1, rat_o3co2, rat_o3co2_1, &
365 : selffac, selffrac, indself, forfac, forfrac, indfor, &
366 : minorfrac, scaleminor, scaleminorn2, indminor, &
367 27648 : fracs, taug)
368 :
369 : ! Combine gaseous and aerosol optical depths, if aerosol active
370 : if (iaer == 0) then
371 : do ig = 1, ngptlw
372 : do k = 1, nlay
373 : taut(k,ig) = taug(k,ig)
374 : end do
375 : end do
376 : else if (iaer == 10) then
377 3898368 : do ig = 1, ngptlw
378 243883008 : do k = 1, nlay
379 243855360 : taut(k,ig) = taug(k,ig) + taua(k,ngb(ig))
380 : end do
381 : end do
382 : end if
383 :
384 : ! Call the radiative transfer routine.
385 :
386 : call rtrnmc(nlay, istart, iend, iout, pz, semiss, &
387 : cldfmc, taucmc, planklay, planklev, plankbnd, &
388 : pwvcm, fracs, taut, &
389 : totuflux, totdflux, fnet, htr, &
390 27648 : totuclfl, totdclfl, fnetc, htrc, totufluxs, totdfluxs )
391 :
392 : ! Transfer up and down fluxes and heating rate to output arrays.
393 : ! Vertical indexing goes from bottom to top
394 :
395 1769472 : do k = 0, nlay
396 1741824 : uflx(iplon,k+1) = totuflux(k)
397 1741824 : dflx(iplon,k+1) = totdflux(k)
398 1741824 : uflxc(iplon,k+1) = totuclfl(k)
399 1741824 : dflxc(iplon,k+1) = totdclfl(k)
400 29611008 : uflxs(:,iplon,k+1) = totufluxs(1:nbndlw,k)
401 29638656 : dflxs(:,iplon,k+1) = totdfluxs(1:nbndlw,k)
402 : end do
403 1771776 : do k = 0, nlay-1
404 1714176 : hr(iplon,k+1) = htr(k)
405 1741824 : hrc(iplon,k+1) = htrc(k)
406 : end do
407 :
408 : end do
409 :
410 2304 : end subroutine rrtmg_lw
411 :
412 : !=========================================================================================
413 :
414 27648 : subroutine inatm(iplon, nlay, icld, iaer, &
415 27648 : play, plev, tlay, tlev, tsfc, h2ovmr, &
416 55296 : o3vmr, co2vmr, ch4vmr, o2vmr, n2ovmr, cfc11vmr, cfc12vmr, &
417 55296 : cfc22vmr, ccl4vmr, emis, &
418 27648 : cldfmcl, taucmcl, ciwpmcl, clwpmcl, reicmcl, relqmcl, tauaer, &
419 55296 : pavel, pz, tavel, tz, tbound, semiss, coldry, &
420 27648 : wkl, wbrodl, wx, pwvcm, &
421 27648 : cldfmc, taucmc, ciwpmc, clwpmc, reicmc, dgesmc, relqmc, taua)
422 :
423 : ! Input atmospheric profile from GCM, and prepare it for use in RRTMG_LW.
424 : ! Set other RRTMG_LW input parameters.
425 :
426 : use parrrtm, only: nbndlw, ngptlw, nmol, maxxsec, mxmol
427 : use rrlw_con, only: grav, avogad
428 : use rrlw_wvn, only: ixindx
429 :
430 : ! ----- Input -----
431 : integer, intent(in) :: iplon ! column loop index
432 : integer, intent(in) :: nlay ! Number of model layers
433 : integer, intent(in) :: icld ! clear/cloud and cloud overlap flag
434 : integer, intent(in) :: iaer ! aerosol option flag
435 :
436 : real(kind=r8), intent(in) :: play(:,:) ! Layer pressures (hPa, mb)
437 : ! Dimensions: (ncol,nlay)
438 : real(kind=r8), intent(in) :: plev(:,:) ! Interface pressures (hPa, mb)
439 : ! Dimensions: (ncol,nlay+1)
440 : real(kind=r8), intent(in) :: tlay(:,:) ! Layer temperatures (K)
441 : ! Dimensions: (ncol,nlay)
442 : real(kind=r8), intent(in) :: tlev(:,:) ! Interface temperatures (K)
443 : ! Dimensions: (ncol,nlay+1)
444 : real(kind=r8), intent(in) :: tsfc(:) ! Surface temperature (K)
445 : ! Dimensions: (ncol)
446 : real(kind=r8), intent(in) :: h2ovmr(:,:) ! H2O volume mixing ratio
447 : ! Dimensions: (ncol,nlay)
448 : real(kind=r8), intent(in) :: o3vmr(:,:) ! O3 volume mixing ratio
449 : ! Dimensions: (ncol,nlay)
450 : real(kind=r8), intent(in) :: co2vmr(:,:) ! CO2 volume mixing ratio
451 : ! Dimensions: (ncol,nlay)
452 : real(kind=r8), intent(in) :: ch4vmr(:,:) ! Methane volume mixing ratio
453 : ! Dimensions: (ncol,nlay)
454 : real(kind=r8), intent(in) :: o2vmr(:,:) ! O2 volume mixing ratio
455 : ! Dimensions: (ncol,nlay)
456 : real(kind=r8), intent(in) :: n2ovmr(:,:) ! Nitrous oxide volume mixing ratio
457 : ! Dimensions: (ncol,nlay)
458 : real(kind=r8), intent(in) :: cfc11vmr(:,:) ! CFC11 volume mixing ratio
459 : ! Dimensions: (ncol,nlay)
460 : real(kind=r8), intent(in) :: cfc12vmr(:,:) ! CFC12 volume mixing ratio
461 : ! Dimensions: (ncol,nlay)
462 : real(kind=r8), intent(in) :: cfc22vmr(:,:) ! CFC22 volume mixing ratio
463 : ! Dimensions: (ncol,nlay)
464 : real(kind=r8), intent(in) :: ccl4vmr(:,:) ! CCL4 volume mixing ratio
465 : ! Dimensions: (ncol,nlay)
466 : real(kind=r8), intent(in) :: emis(:,:) ! Surface emissivity
467 : ! Dimensions: (ncol,nbndlw)
468 :
469 : real(kind=r8), intent(in) :: cldfmcl(:,:,:) ! Cloud fraction
470 : ! Dimensions: (ngptlw,ncol,nlay)
471 : real(kind=r8), intent(in) :: ciwpmcl(:,:,:) ! Cloud ice water path (g/m2)
472 : ! Dimensions: (ngptlw,ncol,nlay)
473 : real(kind=r8), intent(in) :: clwpmcl(:,:,:) ! Cloud liquid water path (g/m2)
474 : ! Dimensions: (ngptlw,ncol,nlay)
475 : real(kind=r8), intent(in) :: reicmcl(:,:) ! Cloud ice effective radius (microns)
476 : ! Dimensions: (ncol,nlay)
477 : real(kind=r8), intent(in) :: relqmcl(:,:) ! Cloud water drop effective radius (microns)
478 : ! Dimensions: (ncol,nlay)
479 : real(kind=r8), intent(in) :: taucmcl(:,:,:) ! Cloud optical depth
480 : ! Dimensions: (ngptlw,ncol,nlay)
481 : real(kind=r8), intent(in) :: tauaer(:,:,:) ! Aerosol optical depth
482 : ! Dimensions: (ncol,nlay,nbndlw)
483 :
484 : ! ----- Output -----
485 : ! Atmosphere
486 : real(kind=r8), intent(out) :: pavel(nlay) ! layer pressures (mb)
487 : ! Dimensions: (nlay)
488 : real(kind=r8), intent(out) :: tavel(nlay) ! layer temperatures (K)
489 : ! Dimensions: (nlay)
490 : real(kind=r8), intent(out) :: pz(0:nlay) ! level (interface) pressures (hPa, mb)
491 : ! Dimensions: (0:nlay)
492 : real(kind=r8), intent(out) :: tz(0:nlay) ! level (interface) temperatures (K)
493 : ! Dimensions: (0:nlay)
494 : real(kind=r8), intent(out) :: tbound ! surface temperature (K)
495 : real(kind=r8), intent(out) :: coldry(nlay) ! dry air column density (mol/cm2)
496 : ! Dimensions: (nlay)
497 : real(kind=r8), intent(out) :: wbrodl(nlay) ! broadening gas column density (mol/cm2)
498 : ! Dimensions: (nlay)
499 : real(kind=r8), intent(out) :: wkl(mxmol,nlay) ! molecular amounts (mol/cm-2)
500 : ! Dimensions: (mxmol,nlay)
501 : real(kind=r8), intent(out) :: wx(maxxsec,nlay) ! cross-section amounts (mol/cm-2)
502 : ! Dimensions: (maxxsec,nlay)
503 : real(kind=r8), intent(out) :: pwvcm ! precipitable water vapor (cm)
504 : real(kind=r8), intent(out) :: semiss(nbndlw) ! lw surface emissivity
505 : ! Dimensions: (nbndlw)
506 :
507 : ! Atmosphere/clouds - cldprop
508 :
509 : real(kind=r8), intent(out) :: cldfmc(ngptlw,nlay) ! cloud fraction [mcica]
510 : ! Dimensions: (ngptlw,nlay)
511 : real(kind=r8), intent(out) :: ciwpmc(ngptlw,nlay) ! cloud ice water path [mcica]
512 : ! Dimensions: (ngptlw,nlay)
513 : real(kind=r8), intent(out) :: clwpmc(ngptlw,nlay) ! cloud liquid water path [mcica]
514 : ! Dimensions: (ngptlw,nlay)
515 : real(kind=r8), intent(out) :: relqmc(nlay) ! liquid particle effective radius (microns)
516 : ! Dimensions: (nlay)
517 : real(kind=r8), intent(out) :: reicmc(nlay) ! ice particle effective radius (microns)
518 : ! Dimensions: (nlay)
519 : real(kind=r8), intent(out) :: dgesmc(nlay) ! ice particle generalized effective size (microns)
520 : ! Dimensions: (nlay)
521 : real(kind=r8), intent(out) :: taucmc(ngptlw,nlay) ! cloud optical depth [mcica]
522 : ! Dimensions: (ngptlw,nlay)
523 : real(kind=r8), intent(out) :: taua(nlay,nbndlw) ! Aerosol optical depth
524 : ! Dimensions: (nlay,nbndlw)
525 :
526 : ! ----- Local -----
527 : real(kind=r8), parameter :: amd = 28.9660_r8 ! Effective molecular weight of dry air (g/mol)
528 : real(kind=r8), parameter :: amw = 18.0160_r8 ! Molecular weight of water vapor (g/mol)
529 :
530 : ! Set molecular weight ratios (for converting mmr to vmr)
531 : ! e.g. h2ovmr = h2ommr * amdw)
532 : real(kind=r8), parameter :: amdw = 1.607793_r8 ! Molecular weight of dry air / water vapor
533 : real(kind=r8), parameter :: amdc = 0.658114_r8 ! Molecular weight of dry air / carbon dioxide
534 : real(kind=r8), parameter :: amdo = 0.603428_r8 ! Molecular weight of dry air / ozone
535 : real(kind=r8), parameter :: amdm = 1.805423_r8 ! Molecular weight of dry air / methane
536 : real(kind=r8), parameter :: amdn = 0.658090_r8 ! Molecular weight of dry air / nitrous oxide
537 : real(kind=r8), parameter :: amdc1 = 0.210852_r8 ! Molecular weight of dry air / CFC11
538 : real(kind=r8), parameter :: amdc2 = 0.239546_r8 ! Molecular weight of dry air / CFC12
539 :
540 : real(kind=r8), parameter :: sbc = 5.67e-08_r8 ! Stefan-Boltzmann constant (W/m2K4)
541 :
542 : integer :: isp, l, ix, n, imol, ib, ig ! Loop indices
543 : real(kind=r8) :: amm, amttl, wvttl, wvsh, summol
544 : integer :: temp
545 : !----------------------------------------------------------------------------
546 :
547 1741824 : reicmc(:) = 0.0_r8
548 1741824 : dgesmc(:) = 0.0_r8
549 1741824 : relqmc(:) = 0.0_r8
550 241726464 : cldfmc(:,:) = 0.0_r8
551 241726464 : taucmc(:,:) = 0.0_r8
552 241726464 : ciwpmc(:,:) = 0.0_r8
553 241726464 : clwpmc(:,:) = 0.0_r8
554 66880512 : wkl(:,:) = 0.0_r8
555 8598528 : wx(:,:) = 0.0_r8
556 27896832 : taua(:,:) = 0.0_r8
557 27648 : amttl = 0.0_r8
558 27648 : wvttl = 0.0_r8
559 :
560 : ! Set surface temperature.
561 27648 : tbound = tsfc(iplon)
562 :
563 : ! Install input GCM arrays into RRTMG_LW arrays for pressure, temperature,
564 : ! and molecular amounts.
565 : ! Pressures are input in mb, or are converted to mb here.
566 : ! Molecular amounts are input in volume mixing ratio, or are converted from
567 : ! mass mixing ratio (or specific humidity for h2o) to volume mixing ratio
568 : ! here. These are then converted to molecular amount (molec/cm2) below.
569 : ! The dry air column COLDRY (in molec/cm2) is calculated from the level
570 : ! pressures, pz (in mb), based on the hydrostatic equation and includes a
571 : ! correction to account for h2o in the layer. The molecular weight of moist
572 : ! air (amm) is calculated for each layer.
573 : ! Note: In RRTMG, layer indexing goes from bottom to top, and coding below
574 : ! assumes GCM input fields are also bottom to top. Input layer indexing
575 : ! from GCM fields should be reversed here if necessary.
576 :
577 27648 : pz(0) = plev(iplon,nlay+1)
578 27648 : tz(0) = tlev(iplon,nlay+1)
579 1741824 : do l = 1, nlay
580 1714176 : pavel(l) = play(iplon,nlay-l+1)
581 1714176 : tavel(l) = tlay(iplon,nlay-l+1)
582 1714176 : pz(l) = plev(iplon,nlay-l+1)
583 1714176 : tz(l) = tlev(iplon,nlay-l+1)
584 1714176 : wkl(1,l) = h2ovmr(iplon,nlay-l+1)
585 1714176 : wkl(2,l) = co2vmr(iplon,nlay-l+1)
586 1714176 : wkl(3,l) = o3vmr(iplon,nlay-l+1)
587 1714176 : wkl(4,l) = n2ovmr(iplon,nlay-l+1)
588 1714176 : wkl(6,l) = ch4vmr(iplon,nlay-l+1)
589 1714176 : wkl(7,l) = o2vmr(iplon,nlay-l+1)
590 :
591 1714176 : amm = (1._r8 - wkl(1,l)) * amd + wkl(1,l) * amw
592 :
593 1714176 : coldry(l) = (pz(l-1)-pz(l)) * 1.e3_r8 * avogad / &
594 1714176 : (1.e2_r8 * grav * amm * (1._r8 + wkl(1,l)))
595 :
596 : ! Set cross section molecule amounts from input; convert to vmr if necessary
597 1714176 : wx(1,l) = ccl4vmr(iplon,nlay-l+1)
598 1714176 : wx(2,l) = cfc11vmr(iplon,nlay-l+1)
599 1714176 : wx(3,l) = cfc12vmr(iplon,nlay-l+1)
600 1741824 : wx(4,l) = cfc22vmr(iplon,nlay-l+1)
601 :
602 : end do
603 :
604 55296 : coldry(nlay) = (pz(nlay-1)) * 1.e3_r8 * avogad / &
605 55296 : (1.e2_r8 * grav * amm * (1._r8 + wkl(1,nlay-1)))
606 :
607 : ! At this point all molecular amounts in wkl and wx are in volume mixing ratio;
608 : ! convert to molec/cm2 based on coldry for use in rrtm. also, compute precipitable
609 : ! water vapor for diffusivity angle adjustments in rtrn and rtrnmr.
610 1741824 : do l = 1, nlay
611 : summol = 0.0_r8
612 11999232 : do imol = 2, nmol
613 11999232 : summol = summol + wkl(imol,l)
614 : end do
615 1714176 : wbrodl(l) = coldry(l) * (1._r8 - summol)
616 13713408 : do imol = 1, nmol
617 13713408 : wkl(imol,l) = coldry(l) * wkl(imol,l)
618 : end do
619 1714176 : amttl = amttl + coldry(l)+wkl(1,l)
620 1714176 : wvttl = wvttl + wkl(1,l)
621 8598528 : do ix = 1,maxxsec
622 8570880 : if (ixindx(ix) .ne. 0) then
623 6856704 : wx(ixindx(ix),l) = coldry(l) * wx(ix,l) * 1.e-20_r8
624 : end if
625 : end do
626 : end do
627 :
628 27648 : wvsh = (amw * wvttl) / (amd * amttl)
629 27648 : pwvcm = wvsh * (1.e3_r8 * pz(0)) / (1.e2_r8 * grav)
630 :
631 : ! Set spectral surface emissivity for each longwave band.
632 470016 : do n = 1, nbndlw
633 470016 : semiss(n) = emis(iplon,n)
634 : ! semiss(n) = 1.0_r8
635 : end do
636 :
637 : ! Transfer aerosol optical properties to RRTM variable;
638 : ! modify to reverse layer indexing here if necessary.
639 27648 : if (iaer >= 1) then
640 1714176 : do l = 1, nlay-1
641 28698624 : do ib = 1, nbndlw
642 28670976 : taua(l,ib) = tauaer(iplon,nlay-l,ib)
643 : end do
644 : end do
645 : end if
646 :
647 : ! Transfer cloud fraction and cloud optical properties to RRTM variables,
648 : ! modify to reverse layer indexing here if necessary.
649 :
650 27648 : if (icld >= 1) then
651 :
652 : ! Move incoming GCM cloud arrays to RRTMG cloud arrays.
653 :
654 1714176 : do l = 1, nlay-1
655 237800448 : do ig = 1, ngptlw
656 236113920 : cldfmc(ig,l) = cldfmcl(ig,iplon,nlay-l)
657 236113920 : taucmc(ig,l) = taucmcl(ig,iplon,nlay-l)
658 236113920 : ciwpmc(ig,l) = ciwpmcl(ig,iplon,nlay-l)
659 237800448 : clwpmc(ig,l) = clwpmcl(ig,iplon,nlay-l)
660 : end do
661 1686528 : reicmc(l) = reicmcl(iplon,nlay-l)
662 1714176 : relqmc(l) = relqmcl(iplon,nlay-l)
663 : end do
664 : end if
665 :
666 27648 : end subroutine inatm
667 :
668 : end module rrtmg_lw_rad
|