num_quad.F90 17 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336337338339340341342343344345346347348349350351352353354355356357358359360361362363364365366367368369370371372373374375376377378379380381382383384385386387388389390391392393394395396397398399400401402403404405406407408409410411412413414415416417418419420421422423424425426427428429430431432433434435436437438439440441442443444445446447448449450451452453454455456457458459460461462463464465466467468469470471472473474475476477478479480481482483484485486487488489490491492493494495496497498499500501502503504505506507508509510511512513514515516517518519520521522523524525526527528529530531532533534535536537538539540541542543544545546547548549550551552553554555556557558559560561562563564565566567568569570571572573574575576
  1. !
  2. ! num_quad - numerical quadrature
  3. !
  4. module num_quad
  5. implicit none
  6. ! --- in/out -------------------
  7. private
  8. public :: IntervalQuad_Lin
  9. public :: IntervalQuad_Cos_Lin
  10. public :: IntervalQuad_Const
  11. public :: IntervalSum
  12. ! --- const ---------------------------------------
  13. character(len=*), parameter :: mname = 'num_quad'
  14. contains
  15. !-------------------------------------------------------------------
  16. !
  17. ! Syntax:
  18. !
  19. ! call IntervalQuad_Lin( x, y, a, b, c, ilast )
  20. !
  21. ! Compute an approximation to the integral:
  22. !
  23. ! b
  24. ! int y(x) dx
  25. ! x=a
  26. !
  27. ! Input is formed by two vectors x and y: y contains the function
  28. ! values to the abscissa values contained in x.
  29. ! The result, stored in c is computed using the trapezoid formula
  30. ! throughout the domain of integration.
  31. ! If a<x(1) or b>x(n) then the function is extended to a, resp. b,
  32. ! assuming a constant value of y(1), resp. y(n).
  33. ! The index ilast should be set to 1 in the first call. If subsequent
  34. ! integrals on the same vector are to be computed, it should be kept
  35. ! unaltered between calls. This will speed up the lookup of values
  36. ! if the various integrals are ascending in x.
  37. !
  38. !--------------------------------------------------------------------
  39. subroutine IntervalQuad_Lin( x,y, a,b, c, ilast, status )
  40. use GO , only : gol, goErr
  41. use num_tools, only : interval
  42. ! --- in/out ---------------------
  43. real, intent(in) :: x(:), y(:)
  44. real, intent(in) :: a, b
  45. real, intent(out) :: c
  46. integer, intent(inout) :: ilast
  47. integer, intent(out) :: status
  48. ! --- const ---------------------
  49. character(len=*), parameter :: rname = mname//'/IntervalQuad_Lin'
  50. ! --- local -------------------
  51. integer :: n
  52. real :: ya, yb
  53. integer :: i, i1, i2, iflag
  54. ! --- begin ---------------------
  55. n = size(x)
  56. if ( size(y) /= n ) then
  57. write (gol,'("x and y should have same lengths:")'); call goErr
  58. write (gol,'(" size(x) : ",i6)') size(x); call goErr
  59. write (gol,'(" size(y) : ",i6)') size(y); call goErr
  60. write (gol,'("in ",a)') rname; call goErr; status=1; return
  61. end if
  62. i1=ilast
  63. call interval(x,a,i1,iflag)
  64. if (iflag == -1) then ! a < x(1) case
  65. c=y(1)*(x(1)-a)
  66. else if (iflag == 1) then ! a > x(n) case
  67. c=-y(n)*(a-x(n))
  68. else
  69. ya=y(i1)+(a-x(i1))*(y(i1+1)-y(i1))/(x(i1+1)-x(i1))
  70. c=-0.5*(a-x(i1))*(y(i1)+ya)
  71. end if
  72. i2=i1
  73. call interval(x,b,i2,iflag)
  74. if (iflag == -1) then ! b < x(1) case
  75. c=c-y(1)*(x(1)-b)
  76. else if (iflag == 1) then ! b > x(n) case
  77. c=c+y(n)*(b-x(n))
  78. else
  79. yb=y(i2)+(b-x(i2))*(y(i2+1)-y(i2))/(x(i2+1)-x(i2))
  80. c=c+0.5*(b-x(i2))*(y(i2)+yb)
  81. end if
  82. if (i1 < i2) then
  83. do i=i1,i2-1
  84. c=c+0.5*(x(i+1)-x(i))*(y(i+1)+y(i))
  85. end do
  86. else
  87. do i=i2,i1-1
  88. c=c-0.5*(x(i+1)-x(i))*(y(i+1)+y(i))
  89. end do
  90. end if
  91. ilast=i2
  92. ! ok
  93. status = 0
  94. end subroutine IntervalQuad_Lin
  95. ! ===
  96. !-------------------------------------------------------------------
  97. !
  98. ! Syntax:
  99. !
  100. ! call lquad_cos( x, y, a, b, c, ilast )
  101. !
  102. ! Compute an approximation to the integral:
  103. !
  104. ! b
  105. ! int y(x) cos(x) dx
  106. ! x=a
  107. !
  108. ! Input is formed by two vectors x and y: y contains the function
  109. ! values to the abscissa values contained in x.
  110. ! The result, stored in c is computed using the trapezoid formula
  111. ! throughout the domain of integration.
  112. ! If a<x(1) or b>x(n) then the function is extended to a, resp. b,
  113. ! assuming a constant value of y(1), resp. y(n).
  114. ! The index ilast should be set to 1 in the first call. If subsequent
  115. ! integrals on the same vector are to be computed, it should be kept
  116. ! unaltered between calls. This will speed up the lookup of values
  117. ! if the various integrals are ascending in x.
  118. !
  119. ! Assuming y(x) = c + dx in the gridbox :
  120. !
  121. ! b b
  122. ! int (c +dx) cos(x) dx = [ sin(x) + d( x sin(x) + cos(x) ) ]
  123. ! x=a x=a
  124. !
  125. !--------------------------------------------------------------------
  126. subroutine IntervalQuad_Cos_Lin( x, y, a, b, c, ilast, status )
  127. use GO , only : gol, goErr
  128. use num_tools, only : interval
  129. ! --- in/out --------------------------------
  130. real, intent(in) :: x(:), y(:)
  131. real, intent(in) :: a, b
  132. real, intent(out) :: c
  133. integer, intent(inout) :: ilast
  134. integer, intent(out) :: status
  135. ! --- const ---------------------
  136. character(len=*), parameter :: rname = mname//'/IntervalQuad_Cos_Lin'
  137. ! --- local ----------------------------------
  138. integer :: n
  139. real :: s, co
  140. integer :: i, i1, i2
  141. integer :: iflag
  142. ! --- begin ----------------------------------
  143. n = size(x)
  144. if ( size(y) /= n ) then
  145. write (gol,'("x and y should have same lengths:")'); call goErr
  146. write (gol,'(" size(x) : ",i6)') size(x); call goErr
  147. write (gol,'(" size(y) : ",i6)') size(y); call goErr
  148. write (gol,'("in ",a)') rname; call goErr; status=1; return
  149. end if
  150. i1=ilast
  151. call interval(x,a,i1,iflag)
  152. if (iflag == -1) then ! a < x(1) case
  153. c = y(1)*(sin(x(1))-sin(a))
  154. else if (iflag == 1) then ! a > x(n) case
  155. c = -y(n)*(sin(a)-sin(x(n)))
  156. else
  157. ! int (co + s*x) cos(x) dx = INt (co cos(x) dx + Int s*x cos(x) dx
  158. ! co sin(x) | + s x sin(x) | + s cosx |
  159. s = (y(i1+1)-y(i1))/(x(i1+1)-x(i1))
  160. co = y(i1) - x(i1)*s
  161. c = (co + s*a) *sin(a) - &
  162. (co + s*x(i1))*sin(x(i1)) + &
  163. s*(cos(a)-cos(x(i1)))
  164. c = -c ! negative because outside a-b
  165. end if
  166. i2=i1
  167. call interval(x,b,i2,iflag)
  168. if (iflag == -1) then ! b < x(1) case
  169. c = c - y(1)*(sin(x(1))-sin(b))
  170. else if (iflag == 1) then ! b > x(n) case
  171. c = c + y(n)*(sin(b)-sin(x(n)))
  172. else
  173. s = (y(i2+1)-y(i2))/(x(i2+1)-x(i2))
  174. co = y(i2) - x(i2)*s
  175. c = c + (co + s*b)*sin(b) - &
  176. (co + s*x(i2)) *sin(x(i2)) + &
  177. s*(cos(b)-cos(x(i2)))
  178. end if
  179. if (i1 < i2) then
  180. do i=i1,i2-1
  181. s = (y(i+1)-y(i))/(x(i+1)-x(i))
  182. co = y(i) - x(i)*s
  183. c = c + (co + s*x(i+1))*sin(x(i+1)) - &
  184. (co + s*x(i)) *sin(x(i)) + &
  185. s*(cos(x(i+1))-cos(x(i)))
  186. end do
  187. else
  188. do i=i2,i1-1
  189. s = (y(i+1)-y(i))/(x(i+1)-x(i))
  190. co = y(i) - x(i)*s
  191. c = c - (co + s*x(i+1))*sin(x(i+1)) + &
  192. (co + s*x(i)) *sin(x(i)) - &
  193. s*(cos(x(i+1))-cos(x(i)))
  194. end do
  195. end if
  196. ilast=i2
  197. ! ok
  198. status = 0
  199. end subroutine IntervalQuad_Cos_Lin
  200. !-------------------------------------------------------------------
  201. !
  202. ! NAME
  203. ! Integral_const
  204. !
  205. ! INTERFACE
  206. ! subroutine IntervalQuad_Const( x, y, a, b, c, ilast )
  207. ! real, intent(in) :: x(:), y(:)
  208. ! real, intent(in) :: a, b
  209. ! real, intent(out) :: c
  210. ! integer, intent(inout) :: ilast
  211. !
  212. ! DESCRIPTION
  213. ! Compute integral over [a,b] for function specified by:
  214. ! o interval boundaries x(1),..,x(n),x(n+1)
  215. ! o function values y(1),..,y(n) ; constant in interval
  216. !
  217. ! y(1) y(2) y(3)
  218. ! +-----o-----+
  219. ! |/ / / / / /+-----o----+
  220. ! +----o----+ / / / / / / / / / /|
  221. ! | / / / / / / / / / / / |
  222. ! |/ / / / / / / / / / / /|
  223. ! ----+------|--+-----------+--------|-+---
  224. ! x(1) x(2) x(3) x(4)
  225. ! a b
  226. !
  227. ! The result is stored in c.
  228. ! If a<x(1) or b>x(n+1) then the function is extended to a, resp. b,
  229. ! assuming a constant value of y(1), resp. y(n).
  230. ! The index ilast should be set to 1 in the first call.
  231. ! If subsequent integrals on the same vector are to be computed,
  232. ! it should be kept unaltered between calls. This will speed up the
  233. ! lookup of values if the various integrals are ascending in x.
  234. !
  235. !--------------------------------------------------------------------
  236. subroutine IntervalQuad_Const( x, y, a,b, c, ilast, status )
  237. use GO , only : gol, goErr
  238. use num_tools, only : interval
  239. ! --- in/out ---------------------
  240. real, intent(in) :: x(:), y(:)
  241. real, intent(in) :: a, b
  242. real, intent(out) :: c
  243. integer, intent(inout) :: ilast
  244. integer, intent(out) :: status
  245. ! --- const ---------------------
  246. character(len=*), parameter :: rname = mname//'/IntervalQuad_Const'
  247. ! --- local -------------------
  248. integer :: n
  249. real :: ya, yb
  250. integer :: i, i_a, i_b, iflag
  251. ! --- begin ---------------------
  252. ! check array sizes:
  253. n = size(y)
  254. if ( size(x) /= n+1 ) then
  255. write (gol,'("x should have one element more than y:")'); call goErr
  256. write (gol,'(" size(x) : ",i6)') size(x); call goErr
  257. write (gol,'(" size(y) : ",i6)') size(y); call goErr
  258. write (gol,'("in ",a)') rname; call goErr; status=1; return
  259. end if
  260. ! check interval:
  261. if ( b < a ) then
  262. write (gol,'("found strange interval [a,b] :")'); call goErr
  263. write (gol,'(" a : ",es12.4)') a; call goErr
  264. write (gol,'(" b : ",es12.4)') b; call goErr
  265. write (gol,'("in ",a)') rname; call goErr; status=1; return
  266. end if
  267. ! fill c with contribution of interval including a:
  268. i_a = ilast
  269. call Interval( x, a, i_a, iflag )
  270. select case ( iflag )
  271. case ( -1 )
  272. ! a < x(1) ; extend y to the left:
  273. c = y(1) * ( x(1) - a )
  274. case ( 0 )
  275. ! x(i_a) < a < x(i_a+1)
  276. c = y(i_a) * ( x(i_a+1) - a )
  277. case ( 1 )
  278. ! a > x(n+1) ; extend y to the right
  279. ! negative contribution to integral:
  280. c = - y(n) * ( a - x(n+1) )
  281. case default
  282. write (gol,'("unsupported iflag from call to Interval : ",i6)') iflag
  283. write (gol,'("in ",a)') rname; call goErr; status=1; return
  284. end select
  285. ! add contributions of interval including b:
  286. i_b = i_a
  287. call Interval( x, b, i_b, iflag )
  288. select case ( iflag )
  289. case ( -1 )
  290. ! b < x(1) ; negative contribution
  291. c = c - y(1) * ( x(1) - b )
  292. case ( 0 )
  293. ! x(i_b) < b < x(i_b+1)
  294. c = c + y(i_b) * ( b - x(i_b) )
  295. case ( 1 )
  296. ! b > x(n+1)
  297. c = c + y(n) * ( b - x(n+1) )
  298. case default
  299. write (gol,'("unsupported iflag from call to Interval : ",i6)') iflag
  300. write (gol,'("in ",a)') rname; call goErr; status=1; return
  301. end select
  302. ! add contributions of intermediate intervals
  303. do i = i_a+1, i_b-1
  304. c = c + y(i) * (x(i+1)-x(i))
  305. end do
  306. ! set index of last interval (including b):
  307. ilast = i_b
  308. ! ok
  309. status = 0
  310. end subroutine IntervalQuad_Const
  311. !-------------------------------------------------------------------
  312. !
  313. ! NAME
  314. ! IntervalSum - add contributions of small intervals
  315. !
  316. ! INTERFACE
  317. ! subroutine IntervalSum( x, y, a, b, c, ilast, fac )
  318. ! real, intent(in) :: x(:), y(:)
  319. ! real, intent(in) :: a, b
  320. ! real, intent(out) :: c
  321. ! integer, intent(inout) :: ilast
  322. ! real, intent(out), optional :: fac(:)
  323. !
  324. ! DESCRIPTION
  325. ! Compute sum of values y over invterval [a,b],
  326. ! for function y specified by:
  327. ! o interval boundaries x(1),..,x(n),x(n+1)
  328. ! o function values y(1),..,y(n) ; constant in interval
  329. !
  330. ! y(1) y(2) y(3)
  331. !
  332. ! 100%
  333. ! +-----o-----+ 60%
  334. ! 30% | ^ +-----o----+
  335. ! +----o----+ | ^
  336. ! ^ | |
  337. ! | | |
  338. ! | | |
  339. ! ----+------|--+-----------+--------|-+---
  340. ! x(1) x(2) x(3) x(4)
  341. ! a b
  342. !
  343. ! Each y(i) contributes to the sum for the fraction of [x(i),x(i+1)]
  344. ! covered by [a,b].
  345. ! The result is stored in c.
  346. ! If a<x(1) or b>x(n+1) an error is issued.
  347. !
  348. ! The index ilast should be set to 1 in the first call.
  349. ! If subsequent integrals on the same vector are to be computed,
  350. ! it should be kept unaltered between calls. This will speed up the
  351. ! lookup of values if the various integrals are ascending in x.
  352. !
  353. ! If b < a, the result is negative.
  354. !
  355. ! If present, the array 'fac' is filt with the factors applied to y(:).
  356. !
  357. !--------------------------------------------------------------------
  358. subroutine IntervalSum( x, y, a_in, b_in, c, ilast, status, fac )
  359. use GO , only : gol, goErr
  360. use num_tools, only : interval
  361. ! --- in/out ---------------------
  362. real, intent(in) :: x(:), y(:)
  363. real, intent(in) :: a_in, b_in
  364. real, intent(out) :: c
  365. integer, intent(inout) :: ilast
  366. integer, intent(inout) :: status
  367. real, intent(out), optional :: fac(:)
  368. ! --- const --------------------------
  369. character(len=*), parameter :: rname = mname//'/IntervalSum'
  370. ! --- local -------------------
  371. integer :: n
  372. real :: ya, yb
  373. real :: f
  374. integer :: i, i_a, i_b, iflag
  375. real :: a, b, plusmin
  376. ! --- begin ---------------------
  377. ! check array sizes:
  378. n = size(y)
  379. if ( size(x) /= n+1 ) then
  380. write (gol,'("x should have one element more than y:")'); call goErr
  381. write (gol,'(" size(x) : ",i6)') size(x); call goErr
  382. write (gol,'(" size(y) : ",i6)') size(y); call goErr
  383. write (gol,'("in ",a)') rname; call goErr; status=1; return
  384. end if
  385. ! check size of factor array if present:
  386. if ( present(fac) ) then
  387. if ( size(fac) /= n ) then
  388. write (gol,'("ERROR in ",a)') rname; status=1; return
  389. write (gol,'("fac should have same size as y:")'); call goErr
  390. write (gol,'(" size(y) : ",i6)') size(y); call goErr
  391. write (gol,'(" size(fac) : ",i6)') size(fac); call goErr
  392. write (gol,'("in ",a)') rname; call goErr; status=1; return
  393. end if
  394. end if
  395. ! increasing or decreasing interval ?
  396. if ( a_in <= b_in ) then
  397. a = a_in
  398. b = b_in
  399. plusmin = 1.0
  400. else
  401. a = b_in
  402. b = a_in
  403. plusmin = -1.0
  404. end if
  405. ! init output factors:
  406. if ( present(fac) ) fac = 0.0
  407. ! fill c with contribution of interval including a:
  408. i_a = ilast
  409. call interval(x,a,i_a,iflag)
  410. select case ( iflag )
  411. case ( -1 )
  412. ! a < x(1)
  413. write (gol,'("interval is partly less than x :")')
  414. write (gol,'(" a : ",es12.4)') a; call goErr
  415. write (gol,'(" x(1) : ",es12.4)') x(1); call goErr
  416. write (gol,'("ERROR in ",a)') rname; status=1; return
  417. case ( 0 )
  418. ! x(i_a) < a < x(i_a+1)
  419. f = ( x(i_a+1) - a ) / ( x(i_a+1) - x(i_a) )
  420. c = y(i_a) * f
  421. if ( present(fac) ) fac(i_a) = fac(i_a) + f
  422. case ( 1 )
  423. ! a > x(n+1)
  424. write (gol,'("interval partly exceeds x :")')
  425. write (gol,'(" a : ",es12.4)') a; call goErr
  426. write (gol,'(" x(n+1) : ",es12.4)') x(n+1); call goErr
  427. write (gol,'("ERROR in ",a)') rname; status=1; return
  428. case default
  429. write (gol,'("unsupported iflag from call to Interval : ",i6)') iflag
  430. write (gol,'("in ",a)') rname; call goErr; status=1; return
  431. end select
  432. ! add contributions of interval including b:
  433. i_b = i_a
  434. call interval( x, b, i_b, iflag )
  435. select case ( iflag )
  436. case ( -1 )
  437. ! b < x(1)
  438. write (gol,'("interval is outside x :")')
  439. write (gol,'(" b : ",es12.4)') a; call goErr
  440. write (gol,'(" x(1) : ",es12.4)') x(1); call goErr
  441. write (gol,'("ERROR in ",a)') rname; status=1; return
  442. case ( 0 )
  443. ! x(i_b) < b < x(i_b+1)
  444. if ( i_b > i_a ) then
  445. ! b in other interval; add contrib [x(i_b),b]
  446. f = ( b - x(i_b) ) / ( x(i_b+1) - x(i_b) )
  447. else
  448. ! a and b in same interval; substract contrib [b,x(i_b+1)]
  449. f = - ( x(i_b+1) - b ) / ( x(i_b+1) - x(i_b) )
  450. end if
  451. c = c + y(i_b) * f
  452. if ( present(fac) ) fac(i_b) = fac(i_b) + f
  453. case ( 1 )
  454. ! b > x(n+1)
  455. write (gol,'("interval exceeds x :")')
  456. write (gol,'(" b : ",es12.4)') b; call goErr
  457. write (gol,'(" x(n+1) : ",es12.4)') x(n+1); call goErr
  458. write (gol,'("ERROR in ",a)') rname; status=1; return
  459. case default
  460. write (gol,'("unsupported iflag from call to Interval : ",i6)') iflag
  461. write (gol,'("in ",a)') rname; call goErr; status=1; return
  462. end select
  463. ! add contributions of intermediate intervals
  464. if ( present(fac) ) then
  465. do i = i_a+1, i_b-1
  466. c = c + y(i)
  467. fac(i) = fac(i) + 1.0
  468. end do
  469. else
  470. do i = i_a+1, i_b-1
  471. c = c + y(i)
  472. end do
  473. end if
  474. ! set index of last interval (including b):
  475. ilast = i_b
  476. ! apply factor for b<a
  477. c = plusmin * c
  478. ! ok
  479. status = 0
  480. end subroutine IntervalSum
  481. end module num_quad