BUGSrad
 All Classes Files Functions Variables
two_rt_lw_ocastrndm.f90
Go to the documentation of this file.
1 
2 
3 subroutine two_rt_lw_gsolap &
4  ( ncol , nlm , mbs , mbir, &
5  ib , cldamt , wc_cld, wc_clr, &
6  asym_cld ,asym_clr,tau_cld , tau_clr, &
7  es , bf , fu_all , fu_clr, &
8  fd_all , fd_clr , c_maximal, cf_max,&
9  cf_random)
10 
11 
12  use kinds
13  implicit none
14 
15 !-----------------------------------------------------------------------
16 ! REFERENCES:
17 ! Alternate cloud overlap calculations
18 
19 
20 ! SUBROUTINES CALLED:
21 ! none.
22 
23 ! FUNCTIONS CALLED:
24 ! none.
25 
26 ! INCLUDED COMMONS:
27 ! none.
28 
29 ! ARGUMENT LIST VARIABLES:
30 ! INPUT ARGUMENTS:
31 ! ----------------
32  integer (kind=int_kind), intent(in):: &
33  ncol, & !Length of sub-domain.
34  nlm, & !Number of layers.
35  mbs, & !Number of SW spectral intervals.
36  mbir, & !Number of IR spectral intervals.
37  ib !Index of spectral interval.
38 
39  real (kind=dbl_kind), intent(in), dimension(ncol) :: &
40  c_maximal
41 
42  real (kind=dbl_kind), intent(in), dimension(ncol,mbir):: &
43  es !Spectral surface emissivity (-).
44 
45  real (kind=dbl_kind), intent(in), dimension(ncol,nlm):: &
46  wc_cld, & !Single scattering albedo (-).
47  wc_clr, & !Clear sky single-scattering albedo
48  asym_cld, & !Asymmetry factor (-).
49  asym_clr, & !Clear sky asymmetry factor
50  tau_cld, & !Optical depth (-).
51  tau_clr, & !Clear sky optical depth
52  cldamt, & !Cloud fraction (-).
53  cf_max, & !Cf in maximally overlapped region (-).
54  cf_random !Cf in randomly overlapped region (-).
55 
56  real (kind=dbl_kind), intent(in), dimension(ncol,nlm+1):: &
57  bf !Planck function (W/m^2).
58 
59 ! OUTPUT ARGUMENTS:
60 ! -----------------
61  real (kind=dbl_kind), intent(out), dimension(ncol,nlm+1):: &
62  fd_all, & !All-sky spectral downward flux (W/m^2).
63  fd_clr, & !Clear sky spectral downward flux
64  fu_all, & !All-sky spectral upward flux (W/m^2).
65  fu_clr !Clear sky spectral upward flux
66 
67 ! LOCAL VARIABLES:
68 ! ----------------
69  integer (kind=int_kind) :: &
70  i, & !Horizontal index.
71  l, & !Vertical index.
72  ibms !Index of spectral interval.
73 
74  real (kind=dbl_kind), dimension(nlm):: &
75  rr_clr, & !
76  rr_cld, & !
77  rr_tmp, & !
78  tr_clr, & !
79  tr_cld, & !
80  tr_tmp, & !
81  sigu_clr, & !
82  sigu_cld, & !
83  sigu_tmp, & !
84  sigd_clr, & !
85  sigd_cld, & !
86  sigd_tmp !
87 
88  real (kind=dbl_kind) :: &
89  aa, & !
90  bb, & !
91  beta0, & !
92  cc, & !
93  diffac, & !
94  denom, & !
95  fact, & !
96  eggtau, & !
97  ggtau, & !
98  kappa, & !
99  oms, & !
100  prop, & !
101  r, & !
102  rinf, & !
103  t, & !
104  taus !
105 
106  data diffac /2./
107 
108  real (kind=dbl_kind), dimension(nlm):: &
109  td, & !
110  vu, & !
111  exptau
112 
113  real (kind=dbl_kind), dimension(nlm+1):: &
114  re, & !
115  vd, & !
116  fd_max, & !Maximally overlapped cloudy sky spectral downward flux
117  fd_rndm, & !Broken sky (randomly overlapped cloudy sky) spectral downward flux
118  fu_max, & !Maximally overlapped cloud sky spectral upward flux
119  fu_rndm !Broken sky (randomly overlapped cloudy sky) spectral upward flux
120 
121 
122  real (kind=dbl_kind) :: &
123  cld_ran, & !The random overlap fraction of cloud
124  cmin !Minimum cloud fraction in a column
125 
126  ibms = ib - mbs
127 !Begin loop over all columns
128  do i = 1, ncol
129 
130  !Clear layer properties
131  do l = 1, nlm
132  fact = asym_clr(i,l)*asym_clr(i,l)
133  oms = ((1.-fact)*wc_clr(i,l))/(1.-fact*wc_clr(i,l))
134  taus = (1.-fact*wc_clr(i,l))*tau_clr(i,l)
135 
136  beta0 = (4.+asym_clr(i,l))/(8.*(1.+asym_clr(i,l)))
137  t = diffac*(1.-oms*(1.-beta0)) !-0.25
138  r = diffac*oms*beta0 !-0.25
139  kappa = sqrt(t**2-r**2)
140  rinf = r/(kappa+t)
141  ggtau = kappa*taus
142  eggtau = exp(-ggtau)
143  denom = (1.-rinf**2*eggtau**2)
144  tr_clr(l) = (1.-rinf**2)*eggtau/denom
145  rr_clr(l) = rinf*(1.-eggtau**2)/denom
146 
147  if(taus .lt. .8e-2) then
148  sigu_clr(l) = 0.5*diffac*(bf(i,l)+bf(i,l+1))*taus
149  sigd_clr(l) = sigu_clr(l)
150  else
151  aa = (t+r)*(1.-rr_clr(l))-(1.+rr_clr(l)-tr_clr(l))/taus
152  bb = -(t+r)*tr_clr(l)+(1.+rr_clr(l)-tr_clr(l))/taus
153  cc = diffac*(1.-oms)/kappa**2
154  sigu_clr(l) = cc*(aa*bf(i,l)+bb*bf(i,l+1))
155  sigd_clr(l) = cc*(bb*bf(i,l)+aa*bf(i,l+1))
156  endif
157  enddo
158 
159  !CLEAR SKY adding
160  re(1) = 0.
161  vd(1) = 0.
162  do l = 1, nlm
163  prop = 1. / (1. - re(l)*rr_clr(l))
164  re(l+1) = rr_clr(l) + tr_clr(l)**2*re(l)*prop
165  vd(l+1) = sigd_clr(l) + (tr_clr(l)*vd(l) &
166  + tr_clr(l)*re(l)*sigu_clr(l))*prop
167  vu(l) = (rr_clr(l)*vd(l) + sigu_clr(l))*prop
168  td(l) = prop
169  enddo
170 
171  !CLEAR SKY flux calculation
172  fu_clr(i,nlm+1) = es(i,ibms)*bf(i,nlm+1)
173  fd_clr(i,1) = 0.
174  do l = nlm+1, 2, -1
175  fd_clr(i,l) = re(l)*fu_clr(i,l) + vd(l)
176  fu_clr(i,l-1) = tr_clr(l-1)*fu_clr(i,l)*td(l-1) + vu(l-1)
177  enddo
178  !END CLEAR SKY CALC
179 
180 
181  !Cloudy layer properties
182  do l = 1, nlm
183  fact = asym_cld(i,l)*asym_cld(i,l)
184  oms = ((1.-fact)*wc_cld(i,l))/(1.-fact*wc_cld(i,l))
185  taus = (1.-fact*wc_cld(i,l))*tau_cld(i,l)
186 
187  beta0 = (4.+asym_cld(i,l))/(8.*(1.+asym_cld(i,l)))
188  t = diffac*(1.-oms*(1.-beta0)) !-0.25
189  r = diffac*oms*beta0 !-0.25
190  kappa = sqrt(t**2-r**2)
191  rinf = r/(kappa+t)
192  ggtau = kappa*taus
193  eggtau = exp(-ggtau)
194  denom = (1.-rinf**2*eggtau**2)
195  tr_cld(l) = (1.-rinf**2)*eggtau/denom
196  rr_cld(l) = rinf*(1.-eggtau**2)/denom
197 
198  if(taus .lt. .8e-2) then
199  sigu_cld(l) = 0.5*diffac*(bf(i,l)+bf(i,l+1))*taus
200  sigd_cld(l) = sigu_cld(l)
201  else
202  aa = (t+r)*(1.-rr_cld(l))-(1.+rr_cld(l)-tr_cld(l))/taus
203  bb = -(t+r)*tr_cld(l)+(1.+rr_cld(l)-tr_cld(l))/taus
204  cc = diffac*(1.-oms)/kappa**2
205  sigu_cld(l) = cc*(aa*bf(i,l)+bb*bf(i,l+1))
206  sigd_cld(l) = cc*(bb*bf(i,l)+aa*bf(i,l+1))
207  endif
208  enddo
209 
210 
211  !BEGIN MAXIMALLY OVERLAPPED CALC
212  !For layers where 0. < cf_max < 1., treat as randomly combined
213  do l = 1, nlm
214  rr_tmp(l) = cf_max(i,l)*rr_cld(l) + (1. - cf_max(i,l))*rr_clr(l)
215  tr_tmp(l) = cf_max(i,l)*tr_cld(l) + (1. - cf_max(i,l))*tr_clr(l)
216  sigu_tmp(l) = cf_max(i,l)*sigu_cld(l) + (1. - cf_max(i,l))*sigu_clr(l)
217  sigd_tmp(l) = cf_max(i,l)*sigd_cld(l) + (1. - cf_max(i,l))*sigd_clr(l)
218  enddo
219  !MAXIMALLY OVERLAPPED adding
220  re(1) = 0.
221  vd(1) = 0.
222  do l = 1, nlm
223  prop = 1. / (1. - re(l)*rr_tmp(l))
224  re(l+1) = rr_tmp(l) + tr_tmp(l)**2*re(l)*prop
225  vd(l+1) = sigd_tmp(l) + (tr_tmp(l)*vd(l) &
226  + tr_tmp(l)*re(l)*sigu_tmp(l))*prop
227  vu(l) = (rr_tmp(l)*vd(l) + sigu_tmp(l))*prop
228  td(l) = prop
229  enddo
230  !MAXIMALLY OVERLAPPED flux calculation
231  fu_max(nlm+1) = es(i,ibms)*bf(i,nlm+1)
232  do l = nlm+1, 2, -1
233  fd_max(l) = re(l)*fu_max(l) + vd(l)
234  fu_max(l-1) = tr_tmp(l-1)*fu_max(l)*td(l-1) + vu(l-1)
235  enddo
236  !END MAXIMALLY-OVERLAPPED CALC
237 
238 
239  !BEGIN RANDOMLY-OVERLAPPED CALC
240  do l = 1, nlm
241  rr_tmp(l) = cf_random(i,l)*rr_cld(l) + (1. - cf_random(i,l))*rr_clr(l)
242  tr_tmp(l) = cf_random(i,l)*tr_cld(l) + (1. - cf_random(i,l))*tr_clr(l)
243  sigu_tmp(l) = cf_random(i,l)*sigu_cld(l) + (1. - cf_random(i,l))*sigu_clr(l)
244  sigd_tmp(l) = cf_random(i,l)*sigd_cld(l) + (1. - cf_random(i,l))*sigd_clr(l)
245  enddo
246 
247  !RANDOMLY-OVERLAPPED adding
248  re(1) = 0.
249  vd(1) = 0.
250  do l = 1,nlm
251  prop = 1./(1. - re(l)*rr_tmp(l))
252  re(l+1) = rr_tmp(l) + tr_tmp(l)*tr_tmp(l)*re(l)*prop
253  vd(l+1) = sigd_tmp(l) + (tr_tmp(l)*vd(l) + tr_tmp(l)*re(l)*sigu_tmp(l))*prop
254  vu(l) = (rr_tmp(l)*vd(l) + sigu_tmp(l))*prop
255  td(l) = prop
256  enddo
257 
258  !RANDOMLY-OVERLAPPED flux calculation
259  fu_rndm(nlm+1) = es(i,ibms)*bf(i,nlm+1)
260  do l = nlm+1, 2, -1
261  fd_rndm(l) = re(l)*fu_rndm(l) + vd(l)
262  fu_rndm(l-1) = tr_tmp(l-1)*fu_rndm(l)*td(l-1) + vu(l-1)
263  enddo
264  !END RANDOMLY-OVERLAPPED CALC
265 
266  !ALL-SKY: Begin
267  !This is the combined maximally- and randomly overlapped portions of the sky
268  fd_all(i,1) = 0.
269  fu_all(i,nlm+1) = (1. - c_maximal(i))*fu_rndm(nlm+1) + c_maximal(i)*fu_max(nlm+1)
270  do l = 1, nlm
271  fd_all(i,l+1) = (1. - c_maximal(i))*fd_rndm(l+1) + c_maximal(i)*fd_max(l+1)
272  fu_all(i,l) = (1. - c_maximal(i))*fu_rndm(l) + c_maximal(i)*fu_max(l)
273  enddo
275 
276  enddo !Loop over columns
277  return
278  end subroutine two_rt_lw_gsolap