PageRenderTime 62ms CodeModel.GetById 32ms RepoModel.GetById 0ms app.codeStats 0ms

/fortran/tquadgsq2d.f

https://bitbucket.org/carter/libqd-reference
FORTRAN Legacy | 416 lines | 265 code | 80 blank | 71 comment | 12 complexity | e1276682a91f6940bd58a31ac9a6dc8b MD5 | raw file
  1. module quadglobal
  2. ! This work was supported by the Director, Office of Science, Division
  3. ! of Mathematical, Information, and Computational Sciences of the
  4. ! U.S. Department of Energy under contract number DE-AC02-05CH11231.
  5. ! This program demonstrates the 2D quadrature routine 'quadgsq2d'.
  6. ! David H. Bailey 2007-03-23
  7. ! The function quadgsq2d is suitable to integrate any function that is
  8. ! continuous, infinitely differentiable and integrable on a 2D finite open
  9. ! interval. Singularities or vertical derivatives are permitted on the
  10. ! boundaries. This can also be used for certain integrals on infinite
  11. ! intervals, by making a suitable change of variable.
  12. ! The function(s) to be integrated is(are) defined in external function
  13. ! subprogram(s) -- see the sample function subprograms below. The name(s) of
  14. ! the function subprogram(s) must be included in appropriate type and external
  15. ! statements in the main program.
  16. ! Inputs set in parameter statement below:
  17. ! ndebug Debug level setting. Default = 2.
  18. ! ndigits1 Primary working precision. With QD precision, set to 64.
  19. ! nepsilson1 Log10 of the desired tolerance. Normally set to - ndigits1.
  20. ! nquadl Max number of phases in quadrature routine; adding 1 increases
  21. ! (possibly doubles) the number of accurate digits in the result,
  22. ! but also roughly *quadruples* the run time.
  23. ! nq2 Space parameter for wk and xk arrays in the calling program. By
  24. ! default it is set to 8 * 2^nq1. Increase nq2 if directed by a
  25. ! message in subroutine initqts.
  26. use qdmodule
  27. implicit none
  28. integer ndebug, nquadl, ndigits1, nepsilon1, nwords1, nq1, nq2, nqmx
  29. parameter (ndebug = 2, nquadl = 6, ndigits1 = 64, nepsilon1 = -64, &
  30. nwords1 = 2, nq1 = nquadl, nq2 = 8 * 2**nq1)
  31. type (qd_real) xk(-nq2:nq2), wk(-nq2:nq2)
  32. end module
  33. ! program tquadgsq2d
  34. subroutine f_main
  35. use qdmodule
  36. use quadglobal
  37. implicit none
  38. integer i, n, n1
  39. double precision dplog10q, d1, d2, second, tm0, tm1
  40. type (qd_real) cat, catalan, err, quadgsq2d, fun01, fun02, fun03, fun04, &
  41. t1, t2, t3, t4, x1, x2, y1, y2
  42. external quadgsq2d, catalan, fun01, fun02, fun03, fun04, second
  43. integer*4 old_cw
  44. ! This line must be present in DD and QD main programs.
  45. call f_fpu_fix_start (old_cw)
  46. write (6, 1) ndigits1, nepsilon1, nquadl
  47. 1 format ('Quadgsq2d test'/'Digits =',i6,' Epsilon =',i6,' Quadlevel =',i6)
  48. ! Initialize quadrature tables wk and xk (weights and abscissas).
  49. tm0 = second ()
  50. call initqgs
  51. tm1 = second ()
  52. write (6, 2) tm1 - tm0
  53. 2 format ('Quadrature initialization completed: cpu time =',f12.6)
  54. cat = catalan ()
  55. ! Begin quadrature tests.
  56. write (6, 11)
  57. 11 format (/ &
  58. 'Problem 1: Int_-1^1 Int_-1^1 1/(1+x^2+y^2) dx dy = 4*log(2+sqrt(3))-2*pi/3')
  59. x1 = -1.d0
  60. x2 = 1.d0
  61. y1 = -1.d0
  62. y2 = 1.d0
  63. tm0 = second ()
  64. t1 = quadgsq2d (fun01, x1, x2, y1, y2)
  65. tm1 = second ()
  66. write (6, 3) tm1 - tm0
  67. 3 format ('Quadrature completed: CPU time =',f12.6/'Result =')
  68. call qdwrite (6, t1)
  69. t2 = 4.d0 * log (2.d0 + sqrt (qdreal (3.d0))) - 2.d0 * qdpi () / 3.d0
  70. call decmdq (t2 - t1, d1, n1)
  71. write (6, 4) d1, n1
  72. 4 format ('Actual error =',f10.6,'x10^',i5)
  73. write (6, 12)
  74. 12 format (/&
  75. 'Problem 2: Int_0^pi Int_0^pi log (2-cos(s)-cos(t)) = 4*pi*cat- pi^2*log(2)')
  76. x1 = 0.d0
  77. x2 = qdpi()
  78. y1 = 0.d0
  79. y2 = qdpi()
  80. tm0 = second ()
  81. t1 = quadgsq2d (fun02, x1, x2, y1, y2)
  82. tm1 = second ()
  83. t2 = 4.d0 * qdpi() * cat - qdpi()**2 * log (qdreal (2.d0))
  84. write (6, 3) tm1 - tm0
  85. call qdwrite (6, t1)
  86. call decmdq (t2 - t1, d1, n1)
  87. write (6, 4) d1, n1
  88. write (6, 13)
  89. 13 format (/&
  90. 'Problem 3: Int_0^inf Int_0^inf sqrt(x^2+xy+y^2) * exp(-x-y) = 1 + 3/4*log(3)')
  91. x1 = 0.d0
  92. x2 = 1.d0
  93. y1 = 0.d0
  94. y2 = 1.d0
  95. tm0 = second ()
  96. t1 = quadgsq2d (fun03, x1, x2, y1, y2)
  97. tm1 = second ()
  98. t2 = 1.d0 + 0.75d0 * log (qdreal (3.d0))
  99. write (6, 3) tm1 - tm0
  100. call qdwrite (6, t1)
  101. call decmdq (t2 - t1, d1, n1)
  102. write (6, 4) d1, n1
  103. write (6, 14)
  104. 14 format (/&
  105. 'Problem 4: Int_0^1 Int_0^1 1/(sqrt((1-x)*(1-y))*(x+y)) dx dy = 4*cat')
  106. x1 = 0.d0
  107. x2 = 1.d0
  108. y1 = 0.d0
  109. y2 = 1.d0
  110. tm0 = second ()
  111. t1 = quadgsq2d (fun04, x1, x2, y1, y2)
  112. tm1 = second ()
  113. t2 = 4.d0 * cat
  114. write (6, 3) tm1 - tm0
  115. call qdwrite (6, t1)
  116. call decmdq (t2 - t1, d1, n1)
  117. write (6, 4) d1, n1
  118. call f_fpu_fix_end (old_cw)
  119. stop
  120. end
  121. function fun01 (s, t)
  122. ! fun01 (s,t) = 1/sqrt[1+s^2+t^2]
  123. use qdmodule
  124. implicit none
  125. type (qd_real) fun01, s, t
  126. fun01 = 1.d0 / sqrt (1.d0 + s**2 + t**2)
  127. return
  128. end
  129. function fun02 (s, t)
  130. ! fun02 (s,t) = log (2 - cos(s) - cos(t))
  131. use qdmodule
  132. implicit none
  133. type (qd_real) fun02, s, t, t1
  134. t1 = 2.d0 - cos (s) - cos (t)
  135. if (t1 > 0.d0) then
  136. fun02 = log (2.d0 - cos (s) - cos (t))
  137. else
  138. fun02 = 0.d0
  139. endif
  140. return
  141. end
  142. function fun03 (s, t)
  143. ! fun03 (s,t) = ((1/s-1)^2 + (1/s-1)*(1/t-1) + (1/t-1)^2)
  144. ! / (s^2 * t^2 * exp(1/s + 1/t - 2)
  145. use qdmodule
  146. implicit none
  147. type (qd_real) fun03, s, t, s1, t1, sq
  148. external dplog10q
  149. if (s > 3.d-3 .and. t > 3.d-3) then
  150. s1 = 1.d0 / s - 1.d0
  151. t1 = 1.d0 / t - 1.d0
  152. sq = sqrt (s1**2 + s1 * t1 + t1**2)
  153. fun03 = sq / (s**2 * t**2) * exp (-s1 - t1)
  154. else
  155. fun03 = 0.d0
  156. endif
  157. return
  158. end
  159. function fun04 (s, t)
  160. ! fun04 (s,t) = 1/(sqrt((1-s)*(1-t)) * (s+t))
  161. use qdmodule
  162. implicit none
  163. type (qd_real) fun04, s, t
  164. fun04 = 1.d0 / (sqrt ((1.d0 - s) * (1.d0 - t)) * (s + t))
  165. return
  166. end
  167. subroutine initqgs
  168. ! This subroutine initializes the quadrature arays xk and wk for Gaussian
  169. ! quadrature. It employs a Newton iteration scheme with a dynamic precision
  170. ! level. The argument nq2, which is the space allocated for wk and xk in
  171. ! the calling program, should be at least 8 * 2^nq1 + 100, although a higher
  172. ! value may be required, depending on precision level. Monitor the space
  173. ! figure given in the message below during initialization to be certain.
  174. ! David H Bailey 2002-11-04
  175. use qdmodule
  176. use quadglobal
  177. implicit none
  178. integer i, ierror, is, j, j1, k, n, nwp, nws
  179. double precision pi
  180. parameter (pi = 3.141592653589793238d0)
  181. type (qd_real) eps, r, t1, t2, t3, t4, t5
  182. if (ndebug >= 1) then
  183. write (6, 1)
  184. 1 format ('initqgs: Gaussian quadrature initialization')
  185. endif
  186. eps = 10.d0 ** nepsilon1
  187. wk(0) = 0.d0
  188. xk(0) = 0.d0
  189. n = 3 * 2 ** (nq1 + 1)
  190. write (6, *) 'n, nq2 =', n, nq2
  191. do j = 1, n / 2
  192. ! Compute a double precision estimate of the root.
  193. is = 0
  194. r = cos ((pi * (j - 0.25d0)) / (n + 0.5d0))
  195. ! Compute the j-th root of the n-degree Legendre polynomial using Newton's
  196. ! iteration.
  197. 100 continue
  198. t1 = 1.d0
  199. t2 = 0.d0
  200. do j1 = 1, n
  201. t3 = t2
  202. t2 = t1
  203. t1 = (dble (2 * j1 - 1) * r * t2 - dble (j1 - 1) * t3) / dble (j1)
  204. enddo
  205. t4 = dble (n) * (r * t1 - t2) / (r ** 2 - 1.d0)
  206. t5 = r
  207. r = r - t1 / t4
  208. ! Once convergence is achieved at nwp = 3, then start doubling (almost) the
  209. ! working precision level at each iteration until full precision is reached.
  210. if (abs (r - t5) > eps) goto 100
  211. i = i + 1
  212. if (i > nq2) goto 110
  213. xk(i) = r
  214. xk(-i) = -xk(i)
  215. t4 = dble (n) * (r * t1 - t2) / (r ** 2 - 1.d0)
  216. wk(i) = 2.d0 / ((1.d0 - r ** 2) * t4 ** 2)
  217. wk(-i) = wk(i)
  218. enddo
  219. nqmx = i
  220. if (ndebug >= 2) then
  221. write (6, 2) i
  222. 2 format ('initqgs: Table spaced used =',i8)
  223. endif
  224. goto 130
  225. 110 continue
  226. write (6, 3) nq2
  227. 3 format ('initqgsq: Table space parameter is too small; value =',i8)
  228. stop
  229. 130 continue
  230. return
  231. end
  232. function quadgsq2d (fun, x1, x2, y1, y2)
  233. ! This performs 2-D tanh-sinh quadrature. No attempt is made in this code
  234. ! to estimate errors.
  235. ! David H. Bailey 2007-03-23
  236. use qdmodule
  237. use quadglobal
  238. implicit none
  239. integer i, j, k, n
  240. double precision h
  241. type (qd_real) ax, bx, ay, by, quadgsq2d, fun, s1, t1, t2, &
  242. x1, x2, xx1, y1, y2, yy1
  243. external fun
  244. ax = 0.5d0 * (x2 - x1)
  245. bx = 0.5d0 * (x2 + x1)
  246. ay = 0.5d0 * (y2 - y1)
  247. by = 0.5d0 * (y2 + y1)
  248. if (nqmx == 0) then
  249. write (6, 1)
  250. 1 format ('quadgsq2d: quadrature arrays have not been initialized')
  251. stop
  252. endif
  253. h = 0.5d0 ** nq1
  254. s1 = 0.d0
  255. do k = -nqmx, nqmx
  256. ! write (6, *) k, nqmx
  257. yy1 = ay * xk(k) + by
  258. do j = -nqmx, nqmx
  259. xx1 = ax * xk(j) + bx
  260. t1 = fun (xx1, yy1)
  261. s1 = s1 + wk(j) * wk(k) * t1
  262. enddo
  263. enddo
  264. quadgsq2d = ax * ay * s1
  265. return
  266. end
  267. function catalan ()
  268. use qdmodule
  269. implicit none
  270. integer k
  271. real*8 dk, eps
  272. type (qd_real) catalan, c1, c2, c4, c8, r16, t1, t2, t3
  273. type (qd_real) x1, x2, x3, x4, x5, x6
  274. c1 = 1.d0
  275. c2 = 2.d0
  276. c4 = 4.d0
  277. c8 = 8.d0
  278. r16 = 1.d0 / 16.d0
  279. t1 = 0.d0
  280. t2 = 1.d0
  281. eps = 1.d-64
  282. do k = 0, 10000000
  283. dk = k
  284. t3 = t2 * (c8 / (8.d0 * dk + 1.d0) ** 2 + c8 / (8.d0 * dk + 2.d0) ** 2 &
  285. + c4 / (8.d0 * dk + 3.d0) ** 2 - c2 / (8.d0 * dk + 5.d0) ** 2 &
  286. - c2 / (8.d0 * dk + 6.d0) ** 2 - c1 / (8.d0 * dk + 7.d0) ** 2)
  287. t1 = t1 + t3
  288. t2 = r16 * t2
  289. if (t3 < 1.d-5 * eps) goto 100
  290. enddo
  291. write (6, *) 'catalan: error - contact author'
  292. 100 continue
  293. catalan = 1.d0 / 8.d0 * qdpi() * log (c2) + 1.d0 / 16.d0 * t1
  294. return
  295. end
  296. function dplog10q (a)
  297. ! For input MP value a, this routine returns a DP approximation to log10 (a).
  298. use qdmodule
  299. implicit none
  300. integer ia
  301. double precision da, dplog10q, t1
  302. type (qd_real) a
  303. ! call mpmdc (a%mpr, da, ia)
  304. da = a
  305. ia = 0
  306. if (da .eq. 0.d0) then
  307. dplog10q = -9999.d0
  308. else
  309. dplog10q = log10 (abs (da)) + ia * log10 (2.d0)
  310. endif
  311. 100 continue
  312. return
  313. end
  314. subroutine decmdq (a, b, ib)
  315. ! For input MP value a, this routine returns DP b and integer ib such that
  316. ! a = b * 10^ib, with 1 <= abs (b) < 10 for nonzero a.
  317. use qdmodule
  318. implicit none
  319. integer ia, ib
  320. double precision da, b, t1, xlt
  321. parameter (xlt = 0.3010299956639812d0)
  322. type (qd_real) a
  323. ! call mpmdc (a%mpr, da, ia)
  324. da = a
  325. ia = 0
  326. if (da .ne. 0.d0) then
  327. t1 = xlt * ia + log10 (abs (da))
  328. ib = t1
  329. if (t1 .lt. 0.d0) ib = ib - 1
  330. b = sign (10.d0 ** (t1 - ib), da)
  331. else
  332. b = 0.d0
  333. ib = 0
  334. endif
  335. return
  336. end