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