ETS  \$Id: Doxyfile 2162 2020-02-26 14:16:09Z g2dpc $
 All Classes Files Functions Variables Pages
cubint.f90
Go to the documentation of this file.
1 subroutine cubint ( ntab, xtab, ftab, ia_in, ib_in, result, error )
2  USE itm_types, only : r8
3 
4 ! obtained from
5 ! http://people.sc.fsu.edu/~jburkardt/f_src/intlib/intlib.html
6 
7 !*****************************************************************************80
8 !
9 !! CUBINT approximates an integral using cubic interpolation of data.
10 !
11 ! Discussion:
12 !
13 ! The integral to be approximated is
14 !
15 ! Integral ( XTAB(IB) <= X <= XTAB(IA) ) F(X) DX
16 !
17 ! The routine estimates the error in integration.
18 !
19 ! Modified:
20 !
21 ! 10 February 2006
22 !
23 ! Reference:
24 !
25 ! Philip Davis, Philip Rabinowitz,
26 ! Methods of Numerical Integration,
27 ! Second Edition,
28 ! Dover, 2007,
29 ! ISBN: 0486453391,
30 ! LC: QA299.3.D28.
31 !
32 ! Philip Gill, GF Miller,
33 ! An algorithm for the integration of unequally spaced data,
34 ! The Computer Journal,
35 ! Number 15, Number 1, 1972, pages 80-83.
36 !
37 ! Parameters:
38 !
39 ! Input, integer NTAB, the number of tabulated points.
40 ! NTAB must be at least 4.
41 !
42 ! Input, real (R8) XTAB(NTAB), contains the points at which the
43 ! function was tabulated. XTAB should contain distinct
44 ! values, given in ascending order.
45 !
46 ! Input, real (R8) FTAB(NTAB), contains the tabulated function
47 ! values, FTAB(I) = F(XTAB(I)).
48 !
49 ! Input, integer IA, the entry of XTAB at which integration
50 ! is to begin. IA must be no less than 1 and no greater
51 ! than NTAB.
52 !
53 ! Input, integer IB, the entry of XTAB at which integration
54 ! is to end. IB must be no less than 1 and no greater than
55 ! NTAB.
56 !
57 ! Output, real (R8) RESULT, the approximate value of the
58 ! integral from XTAB(IA) to XTAB(IB) of the function.
59 !
60 ! Output, real (R8) ERROR, an estimate of the error in
61 ! integration.
62 !
63  implicit none
64 
65  integer ntab
66 
67  real (R8) c
68  real (R8) d1
69  real (R8) d2
70  real (R8) d3
71  real (R8) error
72  real (R8) ftab(ntab)
73  real (R8) h1
74  real (R8) h2
75  real (R8) h3
76  real (R8) h4
77  integer i
78  integer ia, ia_in
79  integer ib, ib_in
80  integer ind
81  integer it
82  integer j
83  integer k
84  real (R8) r1
85  real (R8) r2
86  real (R8) r3
87  real (R8) r4
88  real (R8) result
89  real (R8) s
90  real (R8) term
91  real (R8) xtab(ntab)
92 
93  ia = ia_in
94  ib = ib_in
95  result = 0.0_r8
96  error = 0.0_r8
97 
98  if ( ia == ib ) then
99  return
100  end if
101 
102  if ( ntab < 4 ) then
103  write ( *, '(a)' ) ' '
104  write ( *, '(a)' ) 'CUBINT - Fatal error!'
105  write ( *, '(a,i8)' ) ' NTAB must be at least 4, but input NTAB = ', ntab
106  stop
107  end if
108 
109  if ( ia < 1 ) then
110  write ( *, '(a)' ) ' '
111  write ( *, '(a)' ) 'CUBINT - Fatal error!'
112  write ( *, '(a,i8)' ) ' IA must be at least 1, but input IA = ', ia
113  stop
114  end if
115 
116  if ( ntab < ia ) then
117  write ( *, '(a)' ) ' '
118  write ( *, '(a)' ) 'CUBINT - Fatal error!'
119  write ( *, '(a,i8)' ) ' IA must be <= NTAB, but input IA = ', ia
120  stop
121  end if
122 
123  if ( ib < 1 ) then
124  write ( *, '(a)' ) ' '
125  write ( *, '(a)' ) 'CUBINT - Fatal error!'
126  write ( *, '(a,i8)' ) ' IB must be at least 1, but input IB = ', ib
127  stop
128  end if
129 
130  if ( ntab < ib ) then
131  write ( *, '(a)' ) ' '
132  write ( *, '(a)' ) 'CUBINT - Fatal error!'
133  write ( *, '(a,i8)' ) ' IB must be <= NTAB, but input IB = ', ib
134  stop
135  end if
136 !
137 ! Temporarily switch IA and IB, and store minus sign in IND
138 ! so that, while integration is carried out from low X's
139 ! to high ones, the sense of the integral is preserved.
140 !
141  if ( ib < ia ) then
142  ind = -1
143  it = ib
144  ib = ia
145  ia = it
146  else
147  ind = 1
148  end if
149 
150  s = 0.0_r8
151  c = 0.0_r8
152  r4 = 0.0_r8
153  j = ntab-2
154  if ( ia < ntab-1 .or. ntab == 4 ) then
155  j = max( 3, ia )
156  end if
157 
158  k = 4
159  if ( 2 < ib .or. ntab == 4 ) then
160  k = min( ntab, ib + 2 ) - 1
161  end if
162 
163  do i = j, k
164 
165  if ( i <= j ) then
166 
167  h2 = xtab(j-1) - xtab(j-2)
168  d3 = ( ftab(j-1) - ftab(j-2) ) / h2
169  h3 = xtab(j) - xtab(j-1)
170  d1 = ( ftab(j) - ftab(j-1) ) / h3
171  h1 = h2 + h3
172  d2 = ( d1 - d3 ) / h1
173  h4 = xtab(j+1) - xtab(j)
174  r1 = ( ftab(j+1) - ftab(j) ) / h4
175  r2 = ( r1 - d1 ) / ( h4 + h3 )
176  h1 = h1 + h4
177  r3 = (r2-d2) / h1
178 
179  if ( ia <= 1 ) then
180  result = h2 * ( ftab(1) + h2 * ( 0.5_r8 * d3 - h2 &
181  * ( d2 / 6.0_r8 -(h2+h3+h3)*r3/12.0_r8)))
182  s = -h2**3 * (h2*(3.0_r8*h2+5.0_r8*h4)+10.0_r8*h3*h1) / 60.0_r8
183  end if
184 
185  else
186 
187  h4 = xtab(i+1) - xtab(i)
188  r1 = ( ftab(i+1) - ftab(i) ) / h4
189  r4 = h4 + h3
190  r2 = ( r1 - d1 ) / r4
191  r4 = r4 + h2
192  r3 = ( r2 - d2 ) / r4
193  r4 = ( r3 - d3 ) / ( r4 + h1 )
194 
195  end if
196 
197  if ( ia < i .and. i <= ib ) then
198 
199  term = h3 * ( ( ftab(i) + ftab(i-1) ) * 0.5_r8 &
200  -h3 * h3 * ( d2 + r2 + ( h2 - h4 ) * r3 ) / 12.0_r8 )
201  result = result + term
202  c = h3**3 * ( 2.0_r8 * h3 * h3 &
203  + 5.0_r8 * ( h3 * ( h4 + h2 ) + 2.0_r8 * h2 * h4 ) ) / 120.0_r8
204  error = error + (c+s)*r4
205 
206  if ( i /= j ) then
207  s = c
208  else
209  s = s + c + c
210  end if
211 
212  else
213 
214  error = error + r4 * s
215 
216  end if
217 
218  if ( k <= i ) then
219 
220  if ( ntab <= ib ) then
221  term = h4 * ( ftab(ntab) - h4 * ( 0.5 * r1 &
222  + h4 * ( r2 / 6.0_r8 + ( h3 + h3 + h4 ) * r3 / 12.0_r8 )))
223  result = result + term
224  error = error - h4**3 * r4 * &
225  ( h4 * ( 3.0_r8 * h4 + 5.0_r8 * h2 ) &
226  + 10.0_r8 * h3 * ( h2 + h3 + h4 ) ) / 60.0_r8
227  end if
228 
229  if ( ntab-1 <= ib ) then
230  error = error + s * r4
231  end if
232 
233  else
234 
235  h1 = h2
236  h2 = h3
237  h3 = h4
238  d1 = r1
239  d2 = r2
240  d3 = r3
241  end if
242 
243  end do
244 !
245 ! Restore original values of IA and IB, reverse signs
246 ! of RESULT and ERROR, to account for integration
247 ! that proceeded from high X to low X.
248 !
249  if ( ind /= 1 ) then
250  it = ib
251  ib = ia
252  ia = it
253  result = -result
254  error = -error
255  end if
256 
257  return
258 end
real(r8) function r2(a, x, xr, xs, yr, ys, psi, psir, F_dia)
subroutine cubint(ntab, xtab, ftab, ia_in, ib_in, result, error)
Definition: cubint.f90:1