Line data Source code
1 : ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_sw/src/rrtmg_sw_cldprmc.f90,v $
2 : ! author: $Author: mike $
3 : ! revision: $Revision: 1.4 $
4 : ! created: $Date: 2008/01/03 21:35:36 $
5 :
6 : module rrtmg_sw_cldprmc
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 : use shr_kind_mod, only: r8 => shr_kind_r8
19 :
20 : use parrrsw, only: ngptsw
21 :
22 : implicit none
23 :
24 : !=========================================================================================
25 : contains
26 : !=========================================================================================
27 :
28 38400 : subroutine cldprmc_sw(ncol,nlayers, inflag, iceflag, liqflag, zcldfmc, &
29 38400 : ciwpmc, clwpmc, reicmc, dgesmc, relqmc, &
30 38400 : ztaormc, ztaucmc, zssacmc, zasmcmc, fsfcmc)
31 :
32 : ! Purpose: Compute the cloud optical properties for each cloudy layer
33 : ! and g-point interval for use by the McICA method.
34 : ! Note: Only inflag = 0 and inflag=2/liqflag=1/iceflag=2,3 are available;
35 : ! (Hu & Stamnes, Key, and Fu) are implemented.
36 :
37 : ! ------- Input -------
38 :
39 : integer, intent(in) :: nlayers ! total number of layers
40 : integer, intent(in) :: ncol ! total number of layers
41 : integer, intent(in) :: inflag ! see definitions
42 : integer, intent(in) :: iceflag ! see definitions
43 : integer, intent(in) :: liqflag ! see definitions
44 :
45 : real(kind=r8), intent(in) :: zcldfmc(ncol,nlayers,ngptsw) ! cloud fraction [mcica]
46 : ! Dimensions: (ncol,nlayers,ngptsw)
47 : real(kind=r8), intent(in) :: ciwpmc(ncol,ngptsw,nlayers) ! cloud ice water path [mcica]
48 : ! Dimensions: (ncol,ngptsw,nlayers)
49 : real(kind=r8), intent(in) :: clwpmc(ncol,ngptsw,nlayers) ! cloud liquid water path [mcica]
50 : ! Dimensions: (ncol,ngptsw,nlayers)
51 : real(kind=r8), intent(in) :: relqmc(ncol,nlayers) ! cloud liquid particle effective radius (microns)
52 : ! Dimensions: (ncol,nlayers)
53 : real(kind=r8), intent(in) :: reicmc(ncol,nlayers) ! cloud ice particle effective radius (microns)
54 : ! Dimensions: (ncol,nlayers)
55 : real(kind=r8), intent(in) :: dgesmc(ncol,nlayers) ! cloud ice particle generalized effective size (microns)
56 : ! Dimensions: (ncol,nlayers)
57 : real(kind=r8), intent(in) :: fsfcmc(ncol,ngptsw,nlayers) ! cloud forward scattering fraction
58 : ! Dimensions: (ncol,ngptsw,nlayers)
59 :
60 : ! ------- Output -------
61 :
62 : real(kind=r8), intent(inout) :: ztaucmc(ncol,nlayers,ngptsw) ! cloud optical depth (delta scaled)
63 : ! Dimensions: (ncol,nlayers,ngptsw)
64 : real(kind=r8), intent(inout) :: zssacmc(ncol,nlayers,ngptsw) ! single scattering albedo (delta scaled)
65 : ! Dimensions: (ncol,nlayers,ngptsw)
66 : real(kind=r8), intent(inout) :: zasmcmc(ncol,nlayers,ngptsw) ! asymmetry parameter (delta scaled)
67 : ! Dimensions: (ncol,nlayers,ngptsw)
68 : real(kind=r8), intent(out) :: ztaormc(ncol,nlayers,ngptsw) ! cloud optical depth (non-delta scaled)
69 : ! Dimensions: (ncol,nlayers,ngptsw)
70 :
71 : ! ------- Local -------
72 :
73 : integer :: lay, ig, iplon
74 : real(kind=r8), parameter :: eps = 1.e-06_r8 ! epsilon
75 : real(kind=r8), parameter :: cldmin = 1.e-80_r8 ! minimum value for cloud quantities
76 76800 : real(kind=r8) :: cwp(ncol) ! total cloud water path
77 :
78 76800 : real(kind=r8), dimension(ncol) :: taucldorig_a, taucloud_a, ssacloud_a, ffp, ffp1, ffpssa
79 : !----------------------------------------------------------------------------
80 :
81 : ! Main layer loop
82 1305600 : do lay = 1, nlayers
83 :
84 : ! Main g-point interval loop
85 143232000 : do ig = 1, ngptsw
86 :
87 1165063680 : do iplon=1, ncol
88 :
89 1021870080 : ztaormc(iplon,lay,ig) = ztaucmc(iplon,lay,ig)
90 1021870080 : cwp(iplon) = ciwpmc(iplon,ig,lay) + clwpmc(iplon,ig,lay)
91 1021870080 : if (zcldfmc(iplon,lay,ig) .ge. cldmin .and. &
92 141926400 : (cwp(iplon) .ge. cldmin .or. ztaucmc(iplon,lay,ig) .ge. cldmin)) then
93 :
94 : ! (inflag=0): Cloud optical properties input directly
95 147100183 : if (inflag .eq. 0) then
96 : ! Cloud optical properties already defined in ztaucmc, zssacmc, zasmcmc are unscaled;
97 : ! Apply delta-M scaling here (using Henyey-Greenstein approximation)
98 147100183 : taucldorig_a(iplon) = ztaucmc(iplon,lay,ig)
99 147100183 : ffp(iplon) = fsfcmc(iplon,ig,lay)
100 147100183 : ffp1(iplon) = 1.0_r8 - ffp(iplon)
101 147100183 : ffpssa(iplon) = 1.0_r8 - ffp(iplon) * zssacmc(iplon,lay,ig)
102 147100183 : ssacloud_a(iplon) = ffp1(iplon) * zssacmc(iplon,lay,ig) / ffpssa(iplon)
103 147100183 : taucloud_a(iplon) = ffpssa(iplon) * taucldorig_a(iplon)
104 :
105 147100183 : ztaormc(iplon,lay,ig) = taucldorig_a(iplon)
106 147100183 : zssacmc(iplon,lay,ig) = ssacloud_a(iplon)
107 147100183 : ztaucmc(iplon,lay,ig) = taucloud_a(iplon)
108 147100183 : zasmcmc(iplon,lay,ig) = (zasmcmc(iplon,lay,ig) - ffp(iplon)) / (ffp1(iplon))
109 :
110 : end if
111 : end if
112 :
113 : end do
114 : end do
115 : end do
116 :
117 38400 : end subroutine cldprmc_sw
118 :
119 : end module rrtmg_sw_cldprmc
|