Line data Source code
1 : ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_sw/src/rrtmg_sw_setcoef.f90,v $
2 : ! author: $Author: mike $
3 : ! revision: $Revision: 1.2 $
4 : ! created: $Date: 2007/08/23 20:40:14 $
5 :
6 : module rrtmg_sw_setcoef
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 rrsw_ref, only: pref, preflog, tref
23 :
24 : implicit none
25 :
26 : contains
27 :
28 : !----------------------------------------------------------------------------
29 486000 : subroutine setcoef_sw(nlayers, pavel, tavel, pz, tz, tbound, coldry, wkl, &
30 243000 : laytrop, layswtch, laylow, jp, jt, jt1, &
31 729000 : co2mult, colch4, colco2, colh2o, colmol, coln2o, &
32 729000 : colo2, colo3, fac00, fac01, fac10, fac11, &
33 243000 : selffac, selffrac, indself, forfac, forfrac, indfor)
34 : !----------------------------------------------------------------------------
35 : !
36 : ! Purpose: For a given atmosphere, calculate the indices and
37 : ! fractions related to the pressure and temperature interpolations.
38 :
39 : ! Modifications:
40 : ! Original: J. Delamere, AER, Inc. (version 2.5, 02/04/01)
41 : ! Revised: Rewritten and adapted to ECMWF F90, JJMorcrette 030224
42 : ! Revised: For uniform rrtmg formatting, MJIacono, Jul 2006
43 :
44 : ! ------ Declarations -------
45 :
46 : ! ----- Input -----
47 : integer, intent(in) :: nlayers ! total number of layers
48 :
49 : real(kind=r8), intent(in) :: pavel(:) ! layer pressures (mb)
50 : ! Dimensions: (nlayers)
51 : real(kind=r8), intent(in) :: tavel(:) ! layer temperatures (K)
52 : ! Dimensions: (nlayers)
53 : real(kind=r8), intent(in) :: pz(0:) ! level (interface) pressures (hPa, mb)
54 : ! Dimensions: (0:nlayers)
55 : real(kind=r8), intent(in) :: tz(0:) ! level (interface) temperatures (K)
56 : ! Dimensions: (0:nlayers)
57 : real(kind=r8), intent(in) :: tbound ! surface temperature (K)
58 : real(kind=r8), intent(in) :: coldry(:) ! dry air column density (mol/cm2)
59 : ! Dimensions: (nlayers)
60 : real(kind=r8), intent(in) :: wkl(:,:) ! molecular amounts (mol/cm-2)
61 : ! Dimensions: (mxmol,nlayers)
62 :
63 : ! ----- Output -----
64 : integer, intent(out) :: laytrop ! tropopause layer index
65 : integer, intent(out) :: layswtch !
66 : integer, intent(out) :: laylow !
67 :
68 : integer, intent(out) :: jp(:) !
69 : ! Dimensions: (nlayers)
70 : integer, intent(out) :: jt(:) !
71 : ! Dimensions: (nlayers)
72 : integer, intent(out) :: jt1(:) !
73 : ! Dimensions: (nlayers)
74 :
75 : real(kind=r8), intent(out) :: colh2o(:) ! column amount (h2o)
76 : ! Dimensions: (nlayers)
77 : real(kind=r8), intent(out) :: colco2(:) ! column amount (co2)
78 : ! Dimensions: (nlayers)
79 : real(kind=r8), intent(out) :: colo3(:) ! column amount (o3)
80 : ! Dimensions: (nlayers)
81 : real(kind=r8), intent(out) :: coln2o(:) ! column amount (n2o)
82 : ! Dimensions: (nlayers)
83 : real(kind=r8), intent(out) :: colch4(:) ! column amount (ch4)
84 : ! Dimensions: (nlayers)
85 : real(kind=r8), intent(out) :: colo2(:) ! column amount (o2)
86 : ! Dimensions: (nlayers)
87 : real(kind=r8), intent(out) :: colmol(:) !
88 : ! Dimensions: (nlayers)
89 : real(kind=r8), intent(out) :: co2mult(:) !
90 : ! Dimensions: (nlayers)
91 :
92 : integer, intent(out) :: indself(:)
93 : ! Dimensions: (nlayers)
94 : integer, intent(out) :: indfor(:)
95 : ! Dimensions: (nlayers)
96 : real(kind=r8), intent(out) :: selffac(:)
97 : ! Dimensions: (nlayers)
98 : real(kind=r8), intent(out) :: selffrac(:)
99 : ! Dimensions: (nlayers)
100 : real(kind=r8), intent(out) :: forfac(:)
101 : ! Dimensions: (nlayers)
102 : real(kind=r8), intent(out) :: forfrac(:)
103 : ! Dimensions: (nlayers)
104 :
105 : real(kind=r8), intent(out) :: & !
106 : fac00(:), fac01(:), & ! Dimensions: (nlayers)
107 : fac10(:), fac11(:)
108 :
109 : ! ----- Local -----
110 :
111 : integer :: indbound
112 : integer :: indlev0
113 : integer :: lay
114 : integer :: jp1
115 :
116 : real(kind=r8) :: stpfac
117 : real(kind=r8) :: tbndfrac
118 : real(kind=r8) :: t0frac
119 : real(kind=r8) :: plog
120 : real(kind=r8) :: fp
121 : real(kind=r8) :: ft
122 : real(kind=r8) :: ft1
123 : real(kind=r8) :: water
124 : real(kind=r8) :: scalefac
125 : real(kind=r8) :: factor
126 : real(kind=r8) :: co2reg
127 : real(kind=r8) :: compfp
128 :
129 :
130 : ! Initializations
131 243000 : stpfac = 296._r8/1013._r8
132 :
133 243000 : indbound = tbound - 159._r8
134 243000 : tbndfrac = tbound - int(tbound)
135 243000 : indlev0 = tz(0) - 159._r8
136 243000 : t0frac = tz(0) - int(tz(0))
137 :
138 243000 : laytrop = 0
139 243000 : layswtch = 0
140 243000 : laylow = 0
141 :
142 : ! Begin layer loop
143 22842000 : do lay = 1, nlayers
144 : ! Find the two reference pressures on either side of the
145 : ! layer pressure. Store them in JP and JP1. Store in FP the
146 : ! fraction of the difference (in ln(pressure)) between these
147 : ! two values that the layer pressure lies.
148 :
149 22599000 : plog = log(pavel(lay))
150 22599000 : jp(lay) = int(36._r8 - 5*(plog+0.04_r8))
151 22599000 : if (jp(lay) .lt. 1) then
152 0 : jp(lay) = 1
153 22599000 : elseif (jp(lay) .gt. 58) then
154 243000 : jp(lay) = 58
155 : endif
156 22599000 : jp1 = jp(lay) + 1
157 22599000 : fp = min(3._r8, max(-2._r8, 5._r8 * (preflog(jp(lay)) - plog)))
158 :
159 : ! Determine, for each reference pressure (JP and JP1), which
160 : ! reference temperature (these are different for each
161 : ! reference pressure) is nearest the layer temperature but does
162 : ! not exceed it. Store these indices in JT and JT1, resp.
163 : ! Store in FT (resp. FT1) the fraction of the way between JT
164 : ! (JT1) and the next highest reference temperature that the
165 : ! layer temperature falls.
166 :
167 22599000 : jt(lay) = int(3._r8 + (tavel(lay)-tref(jp(lay)))/15._r8)
168 22599000 : if (jt(lay) .lt. 1) then
169 418395 : jt(lay) = 1
170 22180605 : elseif (jt(lay) .gt. 4) then
171 211827 : jt(lay) = 4
172 : endif
173 22599000 : ft = min(3._r8, max(-2._r8, ((tavel(lay)-tref(jp(lay)))/15._r8) - float(jt(lay)-3)))
174 22599000 : jt1(lay) = int(3._r8 + (tavel(lay)-tref(jp1))/15._r8)
175 22599000 : if (jt1(lay) .lt. 1) then
176 142454 : jt1(lay) = 1
177 22456546 : elseif (jt1(lay) .gt. 4) then
178 293370 : jt1(lay) = 4
179 : endif
180 22599000 : ft1 = min(3._r8, max(-2._r8, ((tavel(lay)-tref(jp1))/15._r8) - float(jt1(lay)-3)))
181 :
182 22599000 : water = wkl(1,lay)/coldry(lay)
183 22599000 : scalefac = pavel(lay) * stpfac / tavel(lay)
184 :
185 : ! If the pressure is less than ~100mb, perform a different
186 : ! set of species interpolations.
187 :
188 22599000 : if (plog .le. 4.56_r8) go to 5300
189 11635265 : laytrop = laytrop + 1
190 11635265 : if (plog .ge. 6.62_r8) laylow = laylow + 1
191 :
192 : ! Set up factors needed to separately include the water vapor
193 : ! foreign-continuum in the calculation of absorption coefficient.
194 :
195 11635265 : forfac(lay) = scalefac / (1.+water)
196 11635265 : factor = (332.0_r8-tavel(lay))/36.0_r8
197 11635265 : indfor(lay) = min(2, max(1, int(factor)))
198 11635265 : forfrac(lay) = min(3._r8, max(-2._r8, factor - float(indfor(lay))))
199 :
200 : ! Set up factors needed to separately include the water vapor
201 : ! self-continuum in the calculation of absorption coefficient.
202 :
203 11635265 : selffac(lay) = water * forfac(lay)
204 11635265 : factor = (tavel(lay)-188.0_r8)/7.2_r8
205 11635265 : indself(lay) = min(9, max(1, int(factor)-7))
206 11635265 : selffrac(lay) = min(3._r8, max(-2._r8, factor - float(indself(lay) + 7)))
207 :
208 : ! Calculate needed column amounts.
209 :
210 11635265 : colh2o(lay) = 1.e-20_r8 * wkl(1,lay)
211 11635265 : colco2(lay) = 1.e-20_r8 * wkl(2,lay)
212 11635265 : colo3(lay) = 1.e-20_r8 * wkl(3,lay)
213 : ! colo3(lay) = 0._r8
214 : ! colo3(lay) = colo3(lay)/1.16_r8
215 11635265 : coln2o(lay) = 1.e-20_r8 * wkl(4,lay)
216 11635265 : colch4(lay) = 1.e-20_r8 * wkl(6,lay)
217 11635265 : colo2(lay) = 1.e-20_r8 * wkl(7,lay)
218 11635265 : colmol(lay) = 1.e-20_r8 * coldry(lay) + colh2o(lay)
219 : ! colco2(lay) = 0._r8
220 : ! colo3(lay) = 0._r8
221 : ! coln2o(lay) = 0._r8
222 : ! colch4(lay) = 0._r8
223 : ! colo2(lay) = 0._r8
224 : ! colmol(lay) = 0._r8
225 11635265 : if (colco2(lay) .eq. 0._r8) colco2(lay) = 1.e-32_r8 * coldry(lay)
226 11635265 : if (coln2o(lay) .eq. 0._r8) coln2o(lay) = 1.e-32_r8 * coldry(lay)
227 11635265 : if (colch4(lay) .eq. 0._r8) colch4(lay) = 1.e-32_r8 * coldry(lay)
228 11635265 : if (colo2(lay) .eq. 0._r8) colo2(lay) = 1.e-32_r8 * coldry(lay)
229 : ! Using E = 1334.2 cm-1.
230 11635265 : co2reg = 3.55e-24_r8 * coldry(lay)
231 11635265 : co2mult(lay)= (colco2(lay) - co2reg) * &
232 11635265 : 272.63_r8*exp(-1919.4_r8/tavel(lay))/(8.7604e-4_r8*tavel(lay))
233 22599000 : goto 5400
234 :
235 : ! Above laytrop.
236 : 5300 continue
237 :
238 : ! Set up factors needed to separately include the water vapor
239 : ! foreign-continuum in the calculation of absorption coefficient.
240 :
241 10963735 : forfac(lay) = scalefac / (1.+water)
242 10963735 : factor = (tavel(lay)-188.0_r8)/36.0_r8
243 10963735 : indfor(lay) = 3
244 10963735 : forfrac(lay) = min(3._r8, max(-2._r8, factor - 1.0_r8))
245 :
246 : ! Calculate needed column amounts.
247 :
248 10963735 : colh2o(lay) = 1.e-20_r8 * wkl(1,lay)
249 10963735 : colco2(lay) = 1.e-20_r8 * wkl(2,lay)
250 10963735 : colo3(lay) = 1.e-20_r8 * wkl(3,lay)
251 10963735 : coln2o(lay) = 1.e-20_r8 * wkl(4,lay)
252 10963735 : colch4(lay) = 1.e-20_r8 * wkl(6,lay)
253 10963735 : colo2(lay) = 1.e-20_r8 * wkl(7,lay)
254 10963735 : colmol(lay) = 1.e-20_r8 * coldry(lay) + colh2o(lay)
255 10963735 : if (colco2(lay) .eq. 0._r8) colco2(lay) = 1.e-32_r8 * coldry(lay)
256 10963735 : if (coln2o(lay) .eq. 0._r8) coln2o(lay) = 1.e-32_r8 * coldry(lay)
257 10963735 : if (colch4(lay) .eq. 0._r8) colch4(lay) = 1.e-32_r8 * coldry(lay)
258 10963735 : if (colo2(lay) .eq. 0._r8) colo2(lay) = 1.e-32_r8 * coldry(lay)
259 10963735 : co2reg = 3.55e-24_r8 * coldry(lay)
260 10963735 : co2mult(lay)= (colco2(lay) - co2reg) * &
261 10963735 : 272.63_r8*exp(-1919.4_r8/tavel(lay))/(8.7604e-4_r8*tavel(lay))
262 :
263 10963735 : selffac(lay) = 0._r8
264 10963735 : selffrac(lay)= 0._r8
265 10963735 : indself(lay) = 0
266 :
267 : 5400 continue
268 :
269 : ! We have now isolated the layer ln pressure and temperature,
270 : ! between two reference pressures and two reference temperatures
271 : ! (for each reference pressure). We multiply the pressure
272 : ! fraction FP with the appropriate temperature fractions to get
273 : ! the factors that will be needed for the interpolation that yields
274 : ! the optical depths (performed in routines TAUGBn for band n).
275 :
276 22599000 : compfp = 1._r8 - fp
277 22599000 : fac10(lay) = compfp * ft
278 22599000 : fac00(lay) = compfp * (1._r8 - ft)
279 22599000 : fac11(lay) = fp * ft1
280 22842000 : fac01(lay) = fp * (1._r8 - ft1)
281 :
282 : ! End layer loop
283 : enddo
284 :
285 243000 : end subroutine setcoef_sw
286 :
287 : !***************************************************************************
288 1536 : subroutine swatmref
289 : !***************************************************************************
290 :
291 : save
292 :
293 : ! These pressures are chosen such that the ln of the first pressure
294 : ! has only a few non-zero digits (i.e. ln(PREF(1)) = 6.96000) and
295 : ! each subsequent ln(pressure) differs from the previous one by 0.2.
296 :
297 : pref(:) = (/ &
298 : 1.05363e+03_r8,8.62642e+02_r8,7.06272e+02_r8,5.78246e+02_r8,4.73428e+02_r8, &
299 : 3.87610e+02_r8,3.17348e+02_r8,2.59823e+02_r8,2.12725e+02_r8,1.74164e+02_r8, &
300 : 1.42594e+02_r8,1.16746e+02_r8,9.55835e+01_r8,7.82571e+01_r8,6.40715e+01_r8, &
301 : 5.24573e+01_r8,4.29484e+01_r8,3.51632e+01_r8,2.87892e+01_r8,2.35706e+01_r8, &
302 : 1.92980e+01_r8,1.57998e+01_r8,1.29358e+01_r8,1.05910e+01_r8,8.67114e+00_r8, &
303 : 7.09933e+00_r8,5.81244e+00_r8,4.75882e+00_r8,3.89619e+00_r8,3.18993e+00_r8, &
304 : 2.61170e+00_r8,2.13828e+00_r8,1.75067e+00_r8,1.43333e+00_r8,1.17351e+00_r8, &
305 : 9.60789e-01_r8,7.86628e-01_r8,6.44036e-01_r8,5.27292e-01_r8,4.31710e-01_r8, &
306 : 3.53455e-01_r8,2.89384e-01_r8,2.36928e-01_r8,1.93980e-01_r8,1.58817e-01_r8, &
307 : 1.30029e-01_r8,1.06458e-01_r8,8.71608e-02_r8,7.13612e-02_r8,5.84256e-02_r8, &
308 : 4.78349e-02_r8,3.91639e-02_r8,3.20647e-02_r8,2.62523e-02_r8,2.14936e-02_r8, &
309 1536 : 1.75975e-02_r8,1.44076e-02_r8,1.17959e-02_r8,9.65769e-03_r8 /)
310 :
311 : preflog(:) = (/ &
312 : 6.9600e+00_r8, 6.7600e+00_r8, 6.5600e+00_r8, 6.3600e+00_r8, 6.1600e+00_r8, &
313 : 5.9600e+00_r8, 5.7600e+00_r8, 5.5600e+00_r8, 5.3600e+00_r8, 5.1600e+00_r8, &
314 : 4.9600e+00_r8, 4.7600e+00_r8, 4.5600e+00_r8, 4.3600e+00_r8, 4.1600e+00_r8, &
315 : 3.9600e+00_r8, 3.7600e+00_r8, 3.5600e+00_r8, 3.3600e+00_r8, 3.1600e+00_r8, &
316 : 2.9600e+00_r8, 2.7600e+00_r8, 2.5600e+00_r8, 2.3600e+00_r8, 2.1600e+00_r8, &
317 : 1.9600e+00_r8, 1.7600e+00_r8, 1.5600e+00_r8, 1.3600e+00_r8, 1.1600e+00_r8, &
318 : 9.6000e-01_r8, 7.6000e-01_r8, 5.6000e-01_r8, 3.6000e-01_r8, 1.6000e-01_r8, &
319 : -4.0000e-02_r8,-2.4000e-01_r8,-4.4000e-01_r8,-6.4000e-01_r8,-8.4000e-01_r8, &
320 : -1.0400e+00_r8,-1.2400e+00_r8,-1.4400e+00_r8,-1.6400e+00_r8,-1.8400e+00_r8, &
321 : -2.0400e+00_r8,-2.2400e+00_r8,-2.4400e+00_r8,-2.6400e+00_r8,-2.8400e+00_r8, &
322 : -3.0400e+00_r8,-3.2400e+00_r8,-3.4400e+00_r8,-3.6400e+00_r8,-3.8400e+00_r8, &
323 1536 : -4.0400e+00_r8,-4.2400e+00_r8,-4.4400e+00_r8,-4.6400e+00_r8 /)
324 :
325 : ! These are the temperatures associated with the respective
326 : ! pressures for the MLS standard atmosphere.
327 :
328 : tref(:) = (/ &
329 : 2.9420e+02_r8, 2.8799e+02_r8, 2.7894e+02_r8, 2.6925e+02_r8, 2.5983e+02_r8, &
330 : 2.5017e+02_r8, 2.4077e+02_r8, 2.3179e+02_r8, 2.2306e+02_r8, 2.1578e+02_r8, &
331 : 2.1570e+02_r8, 2.1570e+02_r8, 2.1570e+02_r8, 2.1706e+02_r8, 2.1858e+02_r8, &
332 : 2.2018e+02_r8, 2.2174e+02_r8, 2.2328e+02_r8, 2.2479e+02_r8, 2.2655e+02_r8, &
333 : 2.2834e+02_r8, 2.3113e+02_r8, 2.3401e+02_r8, 2.3703e+02_r8, 2.4022e+02_r8, &
334 : 2.4371e+02_r8, 2.4726e+02_r8, 2.5085e+02_r8, 2.5457e+02_r8, 2.5832e+02_r8, &
335 : 2.6216e+02_r8, 2.6606e+02_r8, 2.6999e+02_r8, 2.7340e+02_r8, 2.7536e+02_r8, &
336 : 2.7568e+02_r8, 2.7372e+02_r8, 2.7163e+02_r8, 2.6955e+02_r8, 2.6593e+02_r8, &
337 : 2.6211e+02_r8, 2.5828e+02_r8, 2.5360e+02_r8, 2.4854e+02_r8, 2.4348e+02_r8, &
338 : 2.3809e+02_r8, 2.3206e+02_r8, 2.2603e+02_r8, 2.2000e+02_r8, 2.1435e+02_r8, &
339 : 2.0887e+02_r8, 2.0340e+02_r8, 1.9792e+02_r8, 1.9290e+02_r8, 1.8809e+02_r8, &
340 1536 : 1.8329e+02_r8, 1.7849e+02_r8, 1.7394e+02_r8, 1.7212e+02_r8 /)
341 :
342 1536 : end subroutine swatmref
343 :
344 : end module rrtmg_sw_setcoef
345 :
346 :
|