CCPP SciDoc v7.0.0  v7.0.0
Common Community Physics Package Developed at DTC
 
Loading...
Searching...
No Matches
surface_perturbation.F90
1
4
6
10
11 implicit none
12
13 private
14
15 public cdfnor, ppfbet
16
17 contains
18
19! mg, sfc-perts ***
20
21! the routines below are used in the percentile matching algorithm for the
22! albedo and vegetation fraction perturbations
23
27 subroutine cdfnor(z,cdfz)
28 use machine
29
30 implicit none
31 real(kind=kind_phys), intent(out) :: cdfz
32 real(kind=kind_phys),intent(in) :: z
33! local vars
34 integer iflag
35 real(kind=kind_phys) del,x,cdfx,eps
36
37 eps = 1.0e-5
38
39
40 ! definition of passed parameters !
41 ! z = value for which the normal CDF is to be computed
42 ! eps = the absolute accuracy requirment for the CDF
43 ! iflag = error indicator on output 0->no errors, 1->errorflag from
44 ! cdfgam, 2->errorflag from cdfgam
45 ! cdfz = the CDF of the standard normal distribution evaluated at z
46
47 del = 2.0*eps
48 if (z.eq.0.0) then
49 cdfz = 0.5
50 else
51 x = 0.5*z*z
52 call cdfgam(x,0.5,del,iflag, cdfx)
53 if (iflag.ne.0) return
54 if (z.gt.0.0) then
55 cdfz = 0.5+0.5*cdfx
56 else
57 cdfz = 0.5-0.5*cdfx
58 endif
59 endif
60
61 return
62 end
63
65 subroutine cdfgam(x,alpha,eps,iflag,cdfx)
66 use machine
67
68 implicit none
69 real(kind=kind_phys), intent(out) :: cdfx
70 real(kind=kind_phys),intent(in) :: x, alpha, eps
71! local vars
72 integer iflag,i,j,k, imax
73 logical LL
74 real(kind=kind_phys) dx, dgln, p,u,epsx,pdfl, eta, bl, uflo
75 data imax, uflo / 5000, 1.0e-37 /
76
77
78 ! definition of passed parameters !
79 ! x = value for which the CDF is to be computed
80 ! alpha = parameter of gamma function (>0)
81 ! eps = the absolute accuracy requirment for the CDF
82 ! iflag = error indicator on output 0->no errors, 1->either alpha or eps
83 ! is <= oflo, 2->number of terms evaluated in the infinite series exceeds
84 ! imax.
85 ! cdf = the CDF evaluated at x
86
87 cdfx = 0.0
88
89 if (alpha.le.uflo.or.eps.le.uflo) then
90 iflag=1
91 return
92 endif
93 iflag=0
94
95 ! check for special case of x
96 if (x.le.0) return
97
98 dx = x
99 call dgamln(alpha,dgln)
100 pdfl = (alpha-1.0)*log(dx)-dx-dgln
101 if (pdfl.lt.log(uflo)) then
102 if (x.ge.alpha) cdfx = 1.0
103 else
104 p = alpha
105 u = exp(pdfl)
106 ll = .true.
107 if (x.ge.p) then
108 k = int(p)
109 if (p.le.real(k)) k = k-1
110 eta = p - real(k)
111 call dgamln(eta,dgln)
112 bl = (eta-1)*log(dx)-dx-dgln
113 ll = bl.gt.log(eps)
114 endif
115 epsx = eps/x
116 if (ll) then
117 do i=0,imax
118 if (u.le.epsx*(p-x)) return
119 u = x*u/p
120 cdfx = cdfx+u
121 p = p+1.0
122 enddo
123 iflag = 2
124 else
125 do j=1,k
126 p=p-1.0
127 if (u.le.epsx*(x-p)) continue
128 cdfx = cdfx+u
129 u = p*u/x
130 enddo
131 cdfx = 1.0-cdfx
132 endif
133 endif
134 return
135 end subroutine cdfgam
136
138 subroutine dgamln(x,dgamlnout)
139
140 use machine
141 implicit none
142 real(kind=kind_phys), intent(in) :: x
143 real(kind=kind_phys), intent(out) :: dgamlnout
144! local vars
145 integer i, n
146 real(kind=kind_phys) absacc, b1, b2, b3, b4, b5, b6, b7, b8
147 real(kind=kind_phys) c, dx, q, r, xmin, xn
148 data xmin, absacc / 6.894d0, 1.0e-15 /
149 data c / 0.918938533204672741780329736d0 /
150 data b1 / 0.833333333333333333333333333d-1 /
151 data b2 / - 0.277777777777777777777777778d-2 /
152 data b3 / 0.793650793650793650793650794d-3 /
153 data b4 / - 0.595238095238095238095238095d-3 /
154 data b5 / 0.841750841750841750841750842d-3 /
155 data b6 / - 0.191752691752691752691752692d-2 /
156 data b7 / 0.641025641025641025641025641d-2 /
157 data b8 / - 0.295506535947712418300653595d-1 /
158
159 if (x.le.0.0) stop '*** x<=0.0 in function dgamln ***'
160 dx = x
161 n = max(0,int(xmin - dx + 1.0d0) )
162 xn = dx + n
163 r = 1.0d0/xn
164 q = r*r
165 dgamlnout = r*( b1+q*( b2+q*( b3+q*( b4+q*( b5+q*( b6+q*( b7+q*b8 ) ) ) ) ) ) ) +c + (xn-0.5d0)*log(xn)-xn
166
167 if (n.gt.0) then
168 q = 1.0d0
169 do i=0, n-1
170 q = q*(dx+i)
171 enddo
172 dgamlnout = dgamlnout-log(q)
173 endif
174
175 if (dgamlnout + absacc.eq.dgamlnout) then
176 print *,' ********* WARNING FROM FUNCTION DGAMLN *********'
177 print *,' REQUIRED ABSOLUTE ACCURACY NOT ATTAINED FOR X = ',x
178 endif
179 return
180 end subroutine dgamln
181
185 subroutine ppfbet(pr,p,q,iflag,x)
186 use machine
187 implicit none
188 real(kind=kind_phys), intent(in) :: pr, p, q
189 real(kind=kind_phys), intent(out) :: x
190 ! local variables
191 integer iflag, iter, itmax
192 real(kind=kind_phys) tol, a, b, fa, fb, fc, cdf, tol1
193 real(kind=kind_phys) c, d, e, xm, s, u, v, r, eps
194 data itmax, eps / 50, 1.0e-12 /
195
196 ! Compute beta distribution value corresponding to the
197 ! probability and distribution parameters a,b.
198 !
199 ! pr - a probability value in the interval [0,1]
200 ! p - the first parameter of the beta(p,q) distribution
201 ! q - the second parameter of the beta(p,q) distribution
202 ! iflag - erro indicator in output, 0-no errors, 1,2-error flags
203 ! from subroutine cdfbet, 3- pr<0 or pr>1, 4-p<=0 or
204 ! q<=0, 5-tol<1.E-8, 6-the cdfs at the endpoints have
205 ! the same sign and no value of x is defined, 7-maximum
206 ! iterations exceeded and current value of x returned
207
208 tol = 1.0e-5
209
210
211 iflag = 0
212 if (pr.lt.0.0.or.pr.gt.1.) then
213 iflag = 3
214 return
215 endif
216 if(min(p,q).le.0.) then
217 iflag =4
218 return
219 endif
220 if (tol.lt.1.0e-8) then
221 iflag = 5
222 return
223 endif
224 a = 0.
225 b = 1.
226 fa = -pr
227 fb = 1.-pr
228 if (fb*fa.gt.0.0) then
229 iflag = 6
230 return
231 endif
232
233 fc = fb
234 do iter =1,itmax
235 if (fb*fc.gt.0.) then
236 c=a
237 fc=fa
238 d = b-a
239 e=d
240 endif
241 if (abs(fc).lt.abs(fb)) then
242 a=b
243 b=c
244 c=a
245 fa=fb
246 fb=fc
247 fc=fa
248 endif
249
250 tol1 = 2.*eps*abs(b)+0.5*tol
251 xm = 0.5*(c-b)
252 if (abs(xm).le.tol1.or.fb.eq.0.0) then
253 x=b
254 return
255 endif
256 if (abs(e).ge.tol1.and.abs(fa).gt.abs(fb)) then
257 s = fb/fa
258 if (a.eq.c) then
259 u = 2.0*xm*s
260 v = 1.0-s
261 else
262 v = fa/fc
263 r = fb/fc
264 u = s*(2.0*xm*v*(v-r)-(b-a)*(r-1.0))
265 v = (v-1.0)*(r-1.0)*(s-1.0)
266 endif
267 if (u.gt.0.0) v = -v
268 u = abs(u)
269 if (2.0*u.lt.min(3.0*xm*v-abs(tol1*v),abs(e*v))) then
270 e = d
271 d = u/v
272 else
273 d = xm
274 e = d
275 endif
276
277 else
278
279 d=xm
280 e=d
281 endif
282
283 a = b
284 fa = fb
285 if (abs(d).gt.tol1) then
286 b = b+d
287 else
288 b = b+sign(tol1,xm)
289 endif
290 call cdfbet(b,p,q,eps,iflag,cdf)
291 if (iflag.ne.0) return
292 fb = cdf-pr
293 enddo
294 x = b
295
296 return
297 end subroutine ppfbet
298
302 subroutine cdfbet(x,p,q,eps,iflag,cdfx)
303 use machine
304
305 ! Computes the value of the cumulative beta distribution at a
306 ! single point x, given the distribution parameters p,q.
307 !
308 ! x - value at which the CDF is to be computed
309 ! p - first parameter of the beta function
310 ! q - second parameter of the beta function
311 ! eps - desired absolute accuracy
312
313 implicit none
314 real(kind=kind_phys), intent(in) :: x, p, q, eps
315 real(kind=kind_phys), intent(out) :: cdfx
316 ! local vars
317 integer iflag, jmax, j
318 logical LL
319 real(kind=kind_phys) dp, dq, gamln, yxeps, w, uflo
320 real(kind=kind_phys) xy, yx, pq, qp, pdfl, u, r, v
321 real(kind=kind_phys) tmp
322 data jmax, w, uflo / 5000, 20.0, 1.0e-30 /
323
324 cdfx = 0.0
325
326 if (p.le.uflo.or.q.le.uflo.or.eps.le.uflo) then
327 iflag = 1
328 endif
329 iflag = 0
330
331 if (x.le.0.0) return
332 if (x.ge.1.0) then
333 cdfx=1.0
334 else
335 ll = (p+w).ge.(p+q+2.0*w)*x
336 if (ll) then
337 xy = x
338 yx = 1.-xy
339 pq = p
340 qp = q
341 else
342 yx = x
343 xy = 1.-yx
344 qp = p
345 pq = q
346 endif
347
348 call gmln(pq,tmp)
349 dp = (pq-1.)*log(xy)-tmp
350 call gmln(qp,tmp)
351 dq = (qp-1.)*log(yx)-tmp
352 call gmln(pq+qp,tmp)
353 pdfl = tmp+dp+dq
354
355 if (pdfl.ge.log(uflo)) then
356 u = exp(pdfl)*xy/pq
357 r = xy/yx
358 do while (qp.gt.1.)
359 if (u.le.eps*(1.-(pq+qp)*xy/(pq+1.))) then
360 if (.not.ll) cdfx = 1.-cdfx
361 return
362 endif
363 cdfx = cdfx+u
364 pq = pq+1.
365 qp = qp-1.
366 u = qp*r*u/pq
367 enddo
368 v = yx*u
369 yxeps = yx*eps
370 do j = 0, jmax
371 if (v.le.yxeps) then
372 if (.not.ll) cdfx = 1.-cdfx
373 return
374 endif
375 cdfx = cdfx + v
376 pq = pq+1.
377 v = (pq+qp-1.)*xy*v/pq
378 enddo
379 iflag = 2
380 endif
381 if (.not.ll) cdfx = 1.-cdfx
382 endif
383
384 end subroutine cdfbet
385
389 subroutine gmln(x,y)
390 use machine
391 ! Computes the natural logarithm of the gamma distribution. Users
392 ! can set the absolute accuracy and corresponding xmin.
393
394 implicit none
395 real(kind=kind_phys), intent(in) :: x
396 real(kind=kind_phys), intent(out) :: y
397! local vars
398 integer i, n
399 real(kind=kind_phys) absacc, b1, b2, b3, b4, b5, b6, b7, b8
400 real(kind=kind_phys) c, dx, q, r, xmin, xn
401! data xmin, absacc / 6.894d0, 1.0E-15 /
402 data xmin, absacc / 1.357d0, 1.0e-3 /
403 data c / 0.918938533204672741780329736d0 /
404 data b1 / 0.833333333333333333333333333d-1 /
405 data b2 / - 0.277777777777777777777777778d-2 /
406 data b3 / 0.793650793650793650793650794d-3 /
407 data b4 / - 0.595238095238095238095238095d-3 /
408 data b5 / 0.841750841750841750841750842d-3 /
409 data b6 / - 0.191752691752691752691752692d-2 /
410 data b7 / 0.641025641025641025641025641d-2 /
411 data b8 / - 0.295506535947712418300653595d-1 /
412
413 if (x.le.0.0) stop '*** x<=0.0 in function gamln ***'
414 dx = x
415 n = max(0,int(xmin - dx + 1.0d0) )
416 xn = dx + n
417 r = 1.0d0/xn
418 q = r*r
419 y = r*( b1+q*( b2+q*( b3+q*( b4+q*( b5+q*( b6+q*( b7+q*b8 ) &
420 & )) ) ) ) ) +c + (xn-0.5d0)*log(xn)-xn
421
422 if (n.gt.0) then
423 q = 1.0d0
424 do i=0, n-1
425 q = q*(dx+i)
426 enddo
427 y = y-log(q)
428 endif
429
430 if (y + absacc.eq.y) then
431 print *,' ********* WARNING FROM FUNCTION GAMLN *********'
432 print *,' REQUIRED ABSOLUTE ACCURACY NOT ATTAINED FOR X = ',x
433 endif
434 return
435 end subroutine gmln
436
437! *** mg, sfc perts
438end module surface_perturbation
subroutine dgamln(x, dgamlnout)
subroutine gmln(x, y)
This subroutine computes the natural logarithm of the gamma distribution. Users can set the absolute ...
subroutine, public ppfbet(pr, p, q, iflag, x)
This subroutine computes the beta distribution value that matches the percentile from the random patt...
subroutine cdfgam(x, alpha, eps, iflag, cdfx)
subroutine, public cdfnor(z, cdfz)
This subrtouine calculates the CDF of the standard normal distribution evaluated at z.
subroutine cdfbet(x, p, q, eps, iflag, cdfx)
This subroutine computes the value of the cumulative beta distribution at a single point x,...
This module contains routines used in the percentile matching algorithm for the albedo and vegetation...