Line data Source code
1 : ! path: $Source: /storm/rc1/cvsroot/rc/rrtmg_sw/src/rrtmg_sw_vrtqdr.f90,v $
2 : ! author: $Author: mike $
3 : ! revision: $Revision: 1.2 $
4 : ! created: $Date: 2007/08/23 20:40:15 $
5 : !
6 : module rrtmg_sw_vrtqdr
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 parkind, only: jpim, jprb
23 : use parrrsw, only: ngptsw
24 :
25 : implicit none
26 :
27 : contains
28 :
29 : ! --------------------------------------------------------------------------
30 8601600 : subroutine vrtqdr_sw(ncol, klev, kw, &
31 8601600 : pref, prefd, ptra, ptrad, &
32 8601600 : pdbt, prdnd, prup, prupd, ptdbt, &
33 8601600 : pfd, pfu)
34 : ! --------------------------------------------------------------------------
35 :
36 : ! Purpose: This routine performs the vertical quadrature integration
37 : !
38 : ! Interface: *vrtqdr_sw* is called from *spcvrt_sw* and *spcvmc_sw*
39 : !
40 : ! Modifications.
41 : !
42 : ! Original: H. Barker
43 : ! Revision: Integrated with rrtmg_sw, J.-J. Morcrette, ECMWF, Oct 2002
44 : ! Revision: Reformatted for consistency with rrtmg_lw: MJIacono, AER, Jul 2006
45 : !
46 : !-----------------------------------------------------------------------
47 :
48 : ! ------- Declarations -------
49 :
50 : ! Input
51 : integer, intent (in) :: ncol
52 : integer, intent (in) :: klev ! number of model layers
53 : integer, intent (in) :: kw(ncol) ! g-point index
54 :
55 : real(kind=r8), intent(in) :: pref(ncol,klev+1) ! direct beam reflectivity
56 : ! Dimensions: (nlayers+1)
57 : real(kind=r8), intent(in) :: prefd(ncol,klev+1) ! diffuse beam reflectivity
58 : ! Dimensions: (nlayers+1)
59 : real(kind=r8), intent(in) :: ptra(ncol,klev+1) ! direct beam transmissivity
60 : ! Dimensions: (nlayers+1)
61 : real(kind=r8), intent(in) :: ptrad(ncol,klev+1) ! diffuse beam transmissivity
62 : ! Dimensions: (nlayers+1)
63 :
64 : real(kind=r8), intent(in) :: pdbt(ncol,klev+1)
65 : ! Dimensions: (nlayers+1)
66 : real(kind=r8), intent(in) :: ptdbt(ncol,klev+1)
67 : ! Dimensions: (nlayers+1)
68 :
69 : real(kind=r8), intent(inout) :: prdnd(ncol,klev+1)
70 : ! Dimensions: (nlayers+1)
71 : real(kind=r8), intent(inout) :: prup(ncol,klev+1)
72 : ! Dimensions: (nlayers+1)
73 : real(kind=r8), intent(inout) :: prupd(ncol,klev+1)
74 : ! Dimensions: (nlayers+1)
75 :
76 : ! Output
77 : real(kind=r8), intent(out) :: pfd(ncol,klev+1,ngptsw) ! downwelling flux (W/m2)
78 : ! Dimensions: (nlayers+1,ngptsw)
79 : ! unadjusted for earth/sun distance or zenith angle
80 : real(kind=r8), intent(out) :: pfu(ncol,klev+1,ngptsw) ! upwelling flux (W/m2)
81 : ! Dimensions: (nlayers+1,ngptsw)
82 : ! unadjusted for earth/sun distance or zenith angle
83 :
84 : ! Local
85 :
86 : integer :: ikp, ikx, jk
87 : integer :: icol
88 : real(kind=r8) :: zreflect
89 17203200 : real(kind=r8) :: ztdn(ncol,klev+1)
90 :
91 : ! Definitions
92 : !
93 : ! pref(jk) direct reflectance
94 : ! prefd(jk) diffuse reflectance
95 : ! ptra(jk) direct transmittance
96 : ! ptrad(jk) diffuse transmittance
97 : !
98 : ! pdbt(jk) layer mean direct beam transmittance
99 : ! ptdbt(jk) total direct beam transmittance at levels
100 : !
101 : !-----------------------------------------------------------------------------
102 :
103 : ! Link lowest layer with surface
104 :
105 70533120 : do icol=1,ncol
106 61931520 : zreflect = 1._r8 / (1._r8 - prefd(icol,klev+1) * prefd(icol,klev))
107 : prup(icol,klev) = pref(icol,klev) + (ptrad(icol,klev) * &
108 : ((ptra(icol,klev) - pdbt(icol,klev)) * prefd(icol,klev+1) + &
109 61931520 : pdbt(icol,klev) * pref(icol,klev+1))) * zreflect
110 : prupd(icol,klev) = prefd(icol,klev) + ptrad(icol,klev) * ptrad(icol,klev) * &
111 70533120 : prefd(icol,klev+1) * zreflect
112 :
113 : ! Pass from bottom to top
114 : end do
115 283852800 : do jk = 1,klev-1
116 2265661440 : do icol=1,ncol
117 1981808640 : ikp = klev+1-jk
118 1981808640 : ikx = ikp-1
119 1981808640 : zreflect = 1._r8 / (1._r8 -prupd(icol,ikp) * prefd(icol,ikx))
120 : prup(icol,ikx) = pref(icol,ikx) + (ptrad(icol,ikx) * &
121 : ((ptra(icol,ikx) - pdbt(icol,ikx)) * prupd(icol,ikp) + &
122 1981808640 : pdbt(icol,ikx) * prup(icol,ikp))) * zreflect
123 : prupd(icol,ikx) = prefd(icol,ikx) + ptrad(icol,ikx) * ptrad(icol,ikx) * &
124 2257059840 : prupd(icol,ikp) * zreflect
125 : enddo
126 : enddo
127 :
128 : ! Upper boundary conditions
129 70533120 : do icol=1,ncol
130 61931520 : ztdn(icol,1) = 1._r8
131 61931520 : prdnd(icol,1) = 0._r8
132 61931520 : ztdn(icol,2) = ptra(icol,1)
133 70533120 : prdnd(icol,2) = prefd(icol,1)
134 :
135 : ! Pass from top to bottom
136 : end do
137 283852800 : do jk = 2,klev
138 2265661440 : do icol=1,ncol
139 1981808640 : ikp = jk+1
140 1981808640 : zreflect = 1._r8 / (1._r8 - prefd(icol,jk) * prdnd(icol,jk))
141 1981808640 : ztdn(icol,ikp) = ptdbt(icol,jk) * ptra(icol,jk) + &
142 : (ptrad(icol,jk) * ((ztdn(icol,jk) - ptdbt(icol,jk)) + &
143 1981808640 : ptdbt(icol,jk) * pref(icol,jk) * prdnd(icol,jk))) * zreflect
144 : prdnd(icol,ikp) = prefd(icol,jk) + ptrad(icol,jk) * ptrad(icol,jk) * &
145 2257059840 : prdnd(icol,jk) * zreflect
146 : enddo
147 : end do
148 : ! Up and down-welling fluxes at levels
149 :
150 301056000 : do jk = 1,klev+1
151 2406727680 : do icol=1,ncol
152 2105671680 : zreflect = 1._r8 / (1._r8 - prdnd(icol,jk) * prupd(icol,jk))
153 2105671680 : pfu(icol,jk,kw(icol)) = (ptdbt(icol,jk) * prup(icol,jk) + &
154 2105671680 : (ztdn(icol,jk) - ptdbt(icol,jk)) * prupd(icol,jk)) * zreflect
155 : pfd(icol,jk,kw(icol)) = ptdbt(icol,jk) + (ztdn(icol,jk) - ptdbt(icol,jk)+ &
156 2398126080 : ptdbt(icol,jk) * prup(icol,jk) * prdnd(icol,jk)) * zreflect
157 : enddo
158 : end do
159 8601600 : end subroutine vrtqdr_sw
160 :
161 : end module rrtmg_sw_vrtqdr
|