source: palm/trunk/LIB/rrtmg/rrtmg_lw_cldprop.f90

Last change on this file was 1585, checked in by maronga, 10 years ago

Added support for RRTMG radiation code

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