1 | ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_lw/src/rrtmg_lw_cldprmc.f90,v $ |
---|
2 | ! author: $Author: miacono $ |
---|
3 | ! revision: $Revision: 1.9 $ |
---|
4 | ! created: $Date: 2011/04/08 20:25:00 $ |
---|
5 | ! |
---|
6 | module rrtmg_lw_cldprmc |
---|
7 | |
---|
8 | ! -------------------------------------------------------------------------- |
---|
9 | ! | | |
---|
10 | ! | Copyright 2002-2009, 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 parkind, only : im => kind_im, rb => kind_rb |
---|
21 | use parrrtm, only : ngptlw, nbndlw |
---|
22 | use rrlw_cld, only: abscld1, absliq0, absliq1, & |
---|
23 | absice0, absice1, absice2, absice3 |
---|
24 | use rrlw_wvn, only: ngb |
---|
25 | use rrlw_vsn, only: hvrclc, hnamclc |
---|
26 | |
---|
27 | implicit none |
---|
28 | |
---|
29 | contains |
---|
30 | |
---|
31 | ! ------------------------------------------------------------------------------ |
---|
32 | subroutine cldprmc(nlayers, inflag, iceflag, liqflag, cldfmc, & |
---|
33 | ciwpmc, clwpmc, reicmc, relqmc, ncbands, taucmc) |
---|
34 | ! ------------------------------------------------------------------------------ |
---|
35 | |
---|
36 | ! Purpose: Compute the cloud optical depth(s) for each cloudy layer. |
---|
37 | |
---|
38 | ! ------- Input ------- |
---|
39 | |
---|
40 | integer(kind=im), intent(in) :: nlayers ! total number of layers |
---|
41 | integer(kind=im), intent(in) :: inflag ! see definitions |
---|
42 | integer(kind=im), intent(in) :: iceflag ! see definitions |
---|
43 | integer(kind=im), intent(in) :: liqflag ! see definitions |
---|
44 | |
---|
45 | real(kind=rb), intent(in) :: cldfmc(:,:) ! cloud fraction [mcica] |
---|
46 | ! Dimensions: (ngptlw,nlayers) |
---|
47 | real(kind=rb), intent(in) :: ciwpmc(:,:) ! cloud ice water path [mcica] |
---|
48 | ! Dimensions: (ngptlw,nlayers) |
---|
49 | real(kind=rb), intent(in) :: clwpmc(:,:) ! cloud liquid water path [mcica] |
---|
50 | ! Dimensions: (ngptlw,nlayers) |
---|
51 | real(kind=rb), intent(in) :: relqmc(:) ! liquid particle effective radius (microns) |
---|
52 | ! Dimensions: (nlayers) |
---|
53 | real(kind=rb), intent(in) :: reicmc(:) ! ice particle effective radius (microns) |
---|
54 | ! Dimensions: (nlayers) |
---|
55 | ! specific definition of reicmc depends on setting of iceflag: |
---|
56 | ! iceflag = 0: ice effective radius, r_ec, (Ebert and Curry, 1992), |
---|
57 | ! r_ec must be >= 10.0 microns |
---|
58 | ! iceflag = 1: ice effective radius, r_ec, (Ebert and Curry, 1992), |
---|
59 | ! r_ec range is limited to 13.0 to 130.0 microns |
---|
60 | ! iceflag = 2: ice effective radius, r_k, (Key, Streamer Ref. Manual, 1996) |
---|
61 | ! r_k range is limited to 5.0 to 131.0 microns |
---|
62 | ! iceflag = 3: generalized effective size, dge, (Fu, 1996), |
---|
63 | ! dge range is limited to 5.0 to 140.0 microns |
---|
64 | ! [dge = 1.0315 * r_ec] |
---|
65 | |
---|
66 | ! ------- Output ------- |
---|
67 | |
---|
68 | integer(kind=im), intent(out) :: ncbands ! number of cloud spectral bands |
---|
69 | real(kind=rb), intent(inout) :: taucmc(:,:) ! cloud optical depth [mcica] |
---|
70 | ! Dimensions: (ngptlw,nlayers) |
---|
71 | |
---|
72 | ! ------- Local ------- |
---|
73 | |
---|
74 | integer(kind=im) :: lay ! Layer index |
---|
75 | integer(kind=im) :: ib ! spectral band index |
---|
76 | integer(kind=im) :: ig ! g-point interval index |
---|
77 | integer(kind=im) :: index |
---|
78 | integer(kind=im) :: icb(nbndlw) |
---|
79 | |
---|
80 | real(kind=rb) :: abscoice(ngptlw) ! ice absorption coefficients |
---|
81 | real(kind=rb) :: abscoliq(ngptlw) ! liquid absorption coefficients |
---|
82 | real(kind=rb) :: cwp ! cloud water path |
---|
83 | real(kind=rb) :: radice ! cloud ice effective size (microns) |
---|
84 | real(kind=rb) :: factor ! |
---|
85 | real(kind=rb) :: fint ! |
---|
86 | real(kind=rb) :: radliq ! cloud liquid droplet radius (microns) |
---|
87 | real(kind=rb), parameter :: eps = 1.e-6_rb ! epsilon |
---|
88 | real(kind=rb), parameter :: cldmin = 1.e-20_rb ! minimum value for cloud quantities |
---|
89 | |
---|
90 | ! ------- Definitions ------- |
---|
91 | |
---|
92 | ! Explanation of the method for each value of INFLAG. Values of |
---|
93 | ! 0 or 1 for INFLAG do not distingish being liquid and ice clouds. |
---|
94 | ! INFLAG = 2 does distinguish between liquid and ice clouds, and |
---|
95 | ! requires further user input to specify the method to be used to |
---|
96 | ! compute the aborption due to each. |
---|
97 | ! INFLAG = 0: For each cloudy layer, the cloud fraction and (gray) |
---|
98 | ! optical depth are input. |
---|
99 | ! INFLAG = 1: For each cloudy layer, the cloud fraction and cloud |
---|
100 | ! water path (g/m2) are input. The (gray) cloud optical |
---|
101 | ! depth is computed as in CCM2. |
---|
102 | ! INFLAG = 2: For each cloudy layer, the cloud fraction, cloud |
---|
103 | ! water path (g/m2), and cloud ice fraction are input. |
---|
104 | ! ICEFLAG = 0: The ice effective radius (microns) is input and the |
---|
105 | ! optical depths due to ice clouds are computed as in CCM3. |
---|
106 | ! ICEFLAG = 1: The ice effective radius (microns) is input and the |
---|
107 | ! optical depths due to ice clouds are computed as in |
---|
108 | ! Ebert and Curry, JGR, 97, 3831-3836 (1992). The |
---|
109 | ! spectral regions in this work have been matched with |
---|
110 | ! the spectral bands in RRTM to as great an extent |
---|
111 | ! as possible: |
---|
112 | ! E&C 1 IB = 5 RRTM bands 9-16 |
---|
113 | ! E&C 2 IB = 4 RRTM bands 6-8 |
---|
114 | ! E&C 3 IB = 3 RRTM bands 3-5 |
---|
115 | ! E&C 4 IB = 2 RRTM band 2 |
---|
116 | ! E&C 5 IB = 1 RRTM band 1 |
---|
117 | ! ICEFLAG = 2: The ice effective radius (microns) is input and the |
---|
118 | ! optical properties due to ice clouds are computed from |
---|
119 | ! the optical properties stored in the RT code, |
---|
120 | ! STREAMER v3.0 (Reference: Key. J., Streamer |
---|
121 | ! User's Guide, Cooperative Institute for |
---|
122 | ! Meteorological Satellite Studies, 2001, 96 pp.). |
---|
123 | ! Valid range of values for re are between 5.0 and |
---|
124 | ! 131.0 micron. |
---|
125 | ! ICEFLAG = 3: The ice generalized effective size (dge) is input |
---|
126 | ! and the optical properties, are calculated as in |
---|
127 | ! Q. Fu, J. Climate, (1998). Q. Fu provided high resolution |
---|
128 | ! tables which were appropriately averaged for the |
---|
129 | ! bands in RRTM_LW. Linear interpolation is used to |
---|
130 | ! get the coefficients from the stored tables. |
---|
131 | ! Valid range of values for dge are between 5.0 and |
---|
132 | ! 140.0 micron. |
---|
133 | ! LIQFLAG = 0: The optical depths due to water clouds are computed as |
---|
134 | ! in CCM3. |
---|
135 | ! LIQFLAG = 1: The water droplet effective radius (microns) is input |
---|
136 | ! and the optical depths due to water clouds are computed |
---|
137 | ! as in Hu and Stamnes, J., Clim., 6, 728-742, (1993). |
---|
138 | ! The values for absorption coefficients appropriate for |
---|
139 | ! the spectral bands in RRTM have been obtained for a |
---|
140 | ! range of effective radii by an averaging procedure |
---|
141 | ! based on the work of J. Pinto (private communication). |
---|
142 | ! Linear interpolation is used to get the absorption |
---|
143 | ! coefficients for the input effective radius. |
---|
144 | |
---|
145 | data icb /1,2,3,3,3,4,4,4,5, 5, 5, 5, 5, 5, 5, 5/ |
---|
146 | |
---|
147 | hvrclc = '$Revision: 1.9 $' |
---|
148 | |
---|
149 | ncbands = 1 |
---|
150 | |
---|
151 | ! This initialization is done in rrtmg_lw_subcol.F90. |
---|
152 | ! do lay = 1, nlayers |
---|
153 | ! do ig = 1, ngptlw |
---|
154 | ! taucmc(ig,lay) = 0.0_rb |
---|
155 | ! enddo |
---|
156 | ! enddo |
---|
157 | |
---|
158 | ! Main layer loop |
---|
159 | do lay = 1, nlayers |
---|
160 | |
---|
161 | do ig = 1, ngptlw |
---|
162 | cwp = ciwpmc(ig,lay) + clwpmc(ig,lay) |
---|
163 | if (cldfmc(ig,lay) .ge. cldmin .and. & |
---|
164 | (cwp .ge. cldmin .or. taucmc(ig,lay) .ge. cldmin)) then |
---|
165 | |
---|
166 | ! Ice clouds and water clouds combined. |
---|
167 | if (inflag .eq. 0) then |
---|
168 | ! Cloud optical depth already defined in taucmc, return to main program |
---|
169 | return |
---|
170 | |
---|
171 | elseif(inflag .eq. 1) then |
---|
172 | stop 'INFLAG = 1 OPTION NOT AVAILABLE WITH MCICA' |
---|
173 | ! cwp = ciwpmc(ig,lay) + clwpmc(ig,lay) |
---|
174 | ! taucmc(ig,lay) = abscld1 * cwp |
---|
175 | |
---|
176 | ! Separate treatement of ice clouds and water clouds. |
---|
177 | elseif(inflag .eq. 2) then |
---|
178 | radice = reicmc(lay) |
---|
179 | |
---|
180 | ! Calculation of absorption coefficients due to ice clouds. |
---|
181 | if (ciwpmc(ig,lay) .eq. 0.0_rb) then |
---|
182 | abscoice(ig) = 0.0_rb |
---|
183 | |
---|
184 | elseif (iceflag .eq. 0) then |
---|
185 | if (radice .lt. 10.0_rb) stop 'ICE RADIUS TOO SMALL' |
---|
186 | abscoice(ig) = absice0(1) + absice0(2)/radice |
---|
187 | |
---|
188 | elseif (iceflag .eq. 1) then |
---|
189 | if (radice .lt. 13.0_rb .or. radice .gt. 130._rb) stop & |
---|
190 | 'ICE RADIUS OUT OF BOUNDS' |
---|
191 | ncbands = 5 |
---|
192 | ib = icb(ngb(ig)) |
---|
193 | abscoice(ig) = absice1(1,ib) + absice1(2,ib)/radice |
---|
194 | |
---|
195 | ! For iceflag=2 option, ice particle effective radius is limited to 5.0 to 131.0 microns |
---|
196 | |
---|
197 | elseif (iceflag .eq. 2) then |
---|
198 | if (radice .lt. 5.0_rb .or. radice .gt. 131.0_rb) stop 'ICE RADIUS OUT OF BOUNDS' |
---|
199 | ncbands = 16 |
---|
200 | factor = (radice - 2._rb)/3._rb |
---|
201 | index = int(factor) |
---|
202 | if (index .eq. 43) index = 42 |
---|
203 | fint = factor - real(index) |
---|
204 | ib = ngb(ig) |
---|
205 | abscoice(ig) = & |
---|
206 | absice2(index,ib) + fint * & |
---|
207 | (absice2(index+1,ib) - (absice2(index,ib))) |
---|
208 | |
---|
209 | ! For iceflag=3 option, ice particle generalized effective size is limited to 5.0 to 140.0 microns |
---|
210 | |
---|
211 | elseif (iceflag .eq. 3) then |
---|
212 | if (radice .lt. 5.0_rb .or. radice .gt. 140.0_rb) stop 'ICE GENERALIZED EFFECTIVE SIZE OUT OF BOUNDS' |
---|
213 | ncbands = 16 |
---|
214 | factor = (radice - 2._rb)/3._rb |
---|
215 | index = int(factor) |
---|
216 | if (index .eq. 46) index = 45 |
---|
217 | fint = factor - real(index) |
---|
218 | ib = ngb(ig) |
---|
219 | abscoice(ig) = & |
---|
220 | absice3(index,ib) + fint * & |
---|
221 | (absice3(index+1,ib) - (absice3(index,ib))) |
---|
222 | |
---|
223 | endif |
---|
224 | |
---|
225 | ! Calculation of absorption coefficients due to water clouds. |
---|
226 | if (clwpmc(ig,lay) .eq. 0.0_rb) then |
---|
227 | abscoliq(ig) = 0.0_rb |
---|
228 | |
---|
229 | elseif (liqflag .eq. 0) then |
---|
230 | abscoliq(ig) = absliq0 |
---|
231 | |
---|
232 | elseif (liqflag .eq. 1) then |
---|
233 | radliq = relqmc(lay) |
---|
234 | if (radliq .lt. 2.5_rb .or. radliq .gt. 60._rb) stop & |
---|
235 | 'LIQUID EFFECTIVE RADIUS OUT OF BOUNDS' |
---|
236 | index = int(radliq - 1.5_rb) |
---|
237 | if (index .eq. 0) index = 1 |
---|
238 | if (index .eq. 58) index = 57 |
---|
239 | fint = radliq - 1.5_rb - real(index) |
---|
240 | ib = ngb(ig) |
---|
241 | abscoliq(ig) = & |
---|
242 | absliq1(index,ib) + fint * & |
---|
243 | (absliq1(index+1,ib) - (absliq1(index,ib))) |
---|
244 | endif |
---|
245 | |
---|
246 | taucmc(ig,lay) = ciwpmc(ig,lay) * abscoice(ig) + & |
---|
247 | clwpmc(ig,lay) * abscoliq(ig) |
---|
248 | |
---|
249 | endif |
---|
250 | endif |
---|
251 | enddo |
---|
252 | enddo |
---|
253 | |
---|
254 | end subroutine cldprmc |
---|
255 | |
---|
256 | end module rrtmg_lw_cldprmc |
---|