-
Notifications
You must be signed in to change notification settings - Fork 315
/
integrate.c
2590 lines (2253 loc) · 82.4 KB
/
integrate.c
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
479
480
481
482
483
484
485
486
487
488
489
490
491
492
493
494
495
496
497
498
499
500
501
502
503
504
505
506
507
508
509
510
511
512
513
514
515
516
517
518
519
520
521
522
523
524
525
526
527
528
529
530
531
532
533
534
535
536
537
538
539
540
541
542
543
544
545
546
547
548
549
550
551
552
553
554
555
556
557
558
559
560
561
562
563
564
565
566
567
568
569
570
571
572
573
574
575
576
577
578
579
580
581
582
583
584
585
586
587
588
589
590
591
592
593
594
595
596
597
598
599
600
601
602
603
604
605
606
607
608
609
610
611
612
613
614
615
616
617
618
619
620
621
622
623
624
625
626
627
628
629
630
631
632
633
634
635
636
637
638
639
640
641
642
643
644
645
646
647
648
649
650
651
652
653
654
655
656
657
658
659
660
661
662
663
664
665
666
667
668
669
670
671
672
673
674
675
676
677
678
679
680
681
682
683
684
685
686
687
688
689
690
691
692
693
694
695
696
697
698
699
700
701
702
703
704
705
706
707
708
709
710
711
712
713
714
715
716
717
718
719
720
721
722
723
724
725
726
727
728
729
730
731
732
733
734
735
736
737
738
739
740
741
742
743
744
745
746
747
748
749
750
751
752
753
754
755
756
757
758
759
760
761
762
763
764
765
766
767
768
769
770
771
772
773
774
775
776
777
778
779
780
781
782
783
784
785
786
787
788
789
790
791
792
793
794
795
796
797
798
799
800
801
802
803
804
805
806
807
808
809
810
811
812
813
814
815
816
817
818
819
820
821
822
823
824
825
826
827
828
829
830
831
832
833
834
835
836
837
838
839
840
841
842
843
844
845
846
847
848
849
850
851
852
853
854
855
856
857
858
859
860
861
862
863
864
865
866
867
868
869
870
871
872
873
874
875
876
877
878
879
880
881
882
883
884
885
886
887
888
889
890
891
892
893
894
895
896
897
898
899
900
901
902
903
904
905
906
907
908
909
910
911
912
913
914
915
916
917
918
919
920
921
922
923
924
925
926
927
928
929
930
931
932
933
934
935
936
937
938
939
940
941
942
943
944
945
946
947
948
949
950
951
952
953
954
955
956
957
958
959
960
961
962
963
964
965
966
967
968
969
970
971
972
973
974
975
976
977
978
979
980
981
982
983
984
985
986
987
988
989
990
991
992
993
994
995
996
997
998
999
1000
/*
* R : A Computer Language for Statistical Data Analysis
* Copyright (C) 2001 the R Development Core Team
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 2 of the License, or
* (at your option) any later version.
*
* This program is distributed in the hope that it will be useful,
* but WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the
* GNU General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 59 Temple Place, Suite 330, Boston, MA 02111-1307 USA
*/
#include <Rinternals.h>
#include <R_ext/Boolean.h>
#ifdef Macintosh
#include <fp.h>
#else
#include <math.h>
#endif
#include <Rmath.h>
typedef struct int_struct
{
SEXP f; /* function */
SEXP env; /* where to evaluate the calls */
} int_struct, *IntStruct;
static void Rdqags(IntStruct IS, double *lower, double *upper,
double *epsabs, double *epsrel, double *result,
double *abserr, int *neval, int *ier, int *limit,
int *lenw, int *last, int *iwork, double *work);
static void Rdqagi(IntStruct IS, double *bound, int *inf,
double *epsabs, double *epsrel, double *result,
double *abserr, int *neval, int *ier, int *limit,
int *lenw, int *last, int *iwork, double *work);
static void Rintfn(double *x, int n, IntStruct IS)
{
SEXP args, resultsxp;
int i;
PROTECT(args = allocVector(REALSXP, n));
for(i = 0; i < n; i++) REAL(args)[i] = x[i];
PROTECT(resultsxp = eval(lang2(IS->f , args), IS->env));
if(length(resultsxp) != n)
error("evaluation of function gave a result of wrong length");
for(i = 0; i < n; i++) {
x[i] = REAL(resultsxp)[i];
if(!R_FINITE(x[i]))
error("non-finite function value");
}
UNPROTECT(2);
return;
}
SEXP call_dqags(SEXP args)
{
int_struct is;
SEXP ans, ansnames;
double lower, upper, epsabs, epsrel, result, abserr, *work;
int neval, ier, limit, lenw, last, *iwork;
args = CDR(args);
is.f = CAR(args); args = CDR(args);
is.env = CAR(args); args = CDR(args);
lower = asReal(CAR(args)); args = CDR(args);
upper = asReal(CAR(args)); args = CDR(args);
epsabs = asReal(CAR(args)); args = CDR(args);
epsrel = asReal(CAR(args)); args = CDR(args);
limit = asInteger(CAR(args)); args = CDR(args);
lenw = 4 * limit;
iwork = (int *) R_alloc(limit, sizeof(int));
work = (double *) R_alloc(lenw, sizeof(double));
Rdqags(&is, &lower, &upper, &epsabs, &epsrel, &result,
&abserr, &neval, &ier, &limit, &lenw, &last, iwork, work);
PROTECT(ans = allocVector(VECSXP, 4));
PROTECT(ansnames = allocVector(STRSXP, 4));
SET_STRING_ELT(ansnames, 0, mkChar("value"));
SET_VECTOR_ELT(ans, 0, allocVector(REALSXP, 1));
REAL(VECTOR_ELT(ans, 0))[0] = result;
SET_STRING_ELT(ansnames, 1, mkChar("abs.error"));
SET_VECTOR_ELT(ans, 1, allocVector(REALSXP, 1));
REAL(VECTOR_ELT(ans, 1))[0] = abserr;
SET_STRING_ELT(ansnames, 2, mkChar("subdivisions"));
SET_VECTOR_ELT(ans, 2, allocVector(INTSXP, 1));
INTEGER(VECTOR_ELT(ans, 2))[0] = last;
SET_STRING_ELT(ansnames, 3, mkChar("ierr"));
SET_VECTOR_ELT(ans, 3, allocVector(INTSXP, 1));
INTEGER(VECTOR_ELT(ans, 3))[0] = ier;
setAttrib(ans, R_NamesSymbol, ansnames);
UNPROTECT(2);
return ans;
}
SEXP call_dqagi(SEXP args)
{
int_struct is;
SEXP ans, ansnames;
double bound, epsabs, epsrel, result, abserr, *work;
int inf, neval, ier, limit, lenw, last, *iwork;
args = CDR(args);
is.f = CAR(args); args = CDR(args);
is.env = CAR(args); args = CDR(args);
bound = asReal(CAR(args)); args = CDR(args);
inf = asInteger(CAR(args)); args = CDR(args);
epsabs = asReal(CAR(args)); args = CDR(args);
epsrel = asReal(CAR(args)); args = CDR(args);
limit = asInteger(CAR(args)); args = CDR(args);
lenw = 4 * limit;
iwork = (int *) R_alloc(limit, sizeof(int));
work = (double *) R_alloc(lenw, sizeof(double));
Rdqagi(&is, &bound,&inf,&epsabs,&epsrel,&result,
&abserr,&neval,&ier,&limit,&lenw,&last,iwork,work);
PROTECT(ans = allocVector(VECSXP, 4));
PROTECT(ansnames = allocVector(STRSXP, 4));
SET_STRING_ELT(ansnames, 0, mkChar("value"));
SET_VECTOR_ELT(ans, 0, allocVector(REALSXP, 1));
REAL(VECTOR_ELT(ans, 0))[0] = result;
SET_STRING_ELT(ansnames, 1, mkChar("abs.error"));
SET_VECTOR_ELT(ans, 1, allocVector(REALSXP, 1));
REAL(VECTOR_ELT(ans, 1))[0] = abserr;
SET_STRING_ELT(ansnames, 2, mkChar("subdivisions"));
SET_VECTOR_ELT(ans, 2, allocVector(INTSXP, 1));
INTEGER(VECTOR_ELT(ans, 2))[0] = last;
SET_STRING_ELT(ansnames, 3, mkChar("ierr"));
SET_VECTOR_ELT(ans, 3, allocVector(INTSXP, 1));
INTEGER(VECTOR_ELT(ans, 3))[0] = ier;
setAttrib(ans, R_NamesSymbol, ansnames);
UNPROTECT(2);
return ans;
}
/* f2c-ed translations + modifications of QUADPACK functions from here down */
static void rdqagie(IntStruct IS, double *, int *, double * , double *, int *,
double *, double *, int *,
int *, double *, double *, double *, double *,
int *, int *);
static void rdqk15i(IntStruct IS, double *, int *, double * , double *,
double *, double *, double *, double *);
static void rdqelg(int *, double *, double *, double *, double *, int *);
static void rdqpsrt(int *, int *, int *, double *, double *, int *, int *);
static void rdqagse(IntStruct IS, double *, double *,
double *, double *, int *, double *, double *,
int *, int *, double *, double *, double *,
double *, int *, int *);
static void rdqk21(IntStruct IS, double *, double *, double *, double *, double *, double *);
static void rdqelg(int *, double *, double *, double *, double *, int *);
/* Table of constant values */
static double c_b6 = 0.;
static double c_b7 = 1.;
static
void Rdqagi(IntStruct IS, double *bound, int *inf, double *
epsabs, double *epsrel, double *result, double *abserr,
int *neval, int *ier, int *limit, int *lenw, int *
last, int *iwork, double *work)
{
int l1, l2, l3;
/*
***begin prologue dqagi
***date written 800101 (yymmdd)
***revision date 830518 (yymmdd)
***category no. h2a3a1,h2a4a1
***keywords automatic integrator, infinite intervals,
general-purpose, transformation, extrapolation,
globally adaptive
***author piessens,robert,appl. math. & progr. div. - k.u.leuven
de doncker,elise,appl. math. & progr. div. -k.u.leuven
***purpose the routine calculates an approximation result to a given
integral i = integral of f over (bound,+infinity)
or i = integral of f over (-infinity,bound)
or i = integral of f over (-infinity,+infinity)
hopefully satisfying following claim for accuracy
abs(i-result).le.max(epsabs,epsrel*abs(i)).
***description
integration over infinite intervals
standard fortran subroutine
parameters
on entry
f - double precision
function subprogram defining the integrand
function f(x). the actual name for f needs to be
declared e x t e r n a l in the driver program.
bound - double precision
finite bound of integration range
(has no meaning if interval is doubly-infinite)
inf - int
indicating the kind of integration range involved
inf = 1 corresponds to (bound,+infinity),
inf = -1 to (-infinity,bound),
inf = 2 to (-infinity,+infinity).
epsabs - double precision
absolute accuracy requested
epsrel - double precision
relative accuracy requested
if epsabs.le.0
and epsrel.lt.max(50*rel.mach.acc.,0.5d-28),
the routine will end with ier = 6.
on return
result - double precision
approximation to the integral
abserr - double precision
estimate of the modulus of the absolute error,
which should equal or exceed abs(i-result)
neval - int
number of integrand evaluations
ier - int
ier = 0 normal and reliable termination of the
routine. it is assumed that the requested
accuracy has been achieved.
- ier.gt.0 abnormal termination of the routine. the
estimates for result and error are less
reliable. it is assumed that the requested
accuracy has not been achieved.
error messages
ier = 1 maximum number of subdivisions allowed
has been achieved. one can allow more
subdivisions by increasing the value of
limit (and taking the according dimension
adjustments into account). however, if
this yields no improvement it is advised
to analyze the integrand in order to
determine the integration difficulties. if
the position of a local difficulty can be
determined (e.g. singularity,
discontinuity within the interval) one
will probably gain from splitting up the
interval at this point and calling the
integrator on the subranges. if possible,
an appropriate special-purpose integrator
should be used, which is designed for
handling the type of difficulty involved.
= 2 the occurrence of roundoff error is
detected, which prevents the requested
tolerance from being achieved.
the error may be under-estimated.
= 3 extremely bad integrand behaviour occurs
at some points of the integration
interval.
= 4 the algorithm does not converge.
roundoff error is detected in the
extrapolation table.
it is assumed that the requested tolerance
cannot be achieved, and that the returned
result is the best which can be obtained.
= 5 the integral is probably divergent, or
slowly convergent. it must be noted that
divergence can occur with any other value
of ier.
= 6 the input is invalid, because
(epsabs.le.0 and
epsrel.lt.max(50*rel.mach.acc.,0.5d-28))
or limit.lt.1 or leniw.lt.limit*4.
result, abserr, neval, last are set to
zero. exept when limit or leniw is
invalid, iwork(1), work(limit*2+1) and
work(limit*3+1) are set to zero, work(1)
is set to a and work(limit+1) to b.
dimensioning parameters
limit - int
dimensioning parameter for iwork
limit determines the maximum number of subintervals
in the partition of the given integration interval
(a,b), limit.ge.1.
if limit.lt.1, the routine will end with ier = 6.
lenw - int
dimensioning parameter for work
lenw must be at least limit*4.
if lenw.lt.limit*4, the routine will end
with ier = 6.
last - int
on return, last equals the number of subintervals
produced in the subdivision process, which
determines the number of significant elements
actually in the work arrays.
work arrays
iwork - int
vector of dimension at least limit, the first
k elements of which contain pointers
to the error estimates over the subintervals,
such that work(limit*3+iwork(1)),... ,
work(limit*3+iwork(k)) form a decreasing
sequence, with k = last if last.le.(limit/2+2), and
k = limit+1-last otherwise
work - double precision
vector of dimension at least lenw
on return
work(1), ..., work(last) contain the left
end points of the subintervals in the
partition of (a,b),
work(limit+1), ..., work(limit+last) contain
the right end points,
work(limit*2+1), ...,work(limit*2+last) contain the
integral approximations over the subintervals,
work(limit*3+1), ..., work(limit*3)
contain the error estimates.
***references (none)
***routines called dqagie
***end prologue dqagi */
/* Parameter adjustments */
--iwork;
--work;
/* Function Body */
*ier = 6;
*neval = 0;
*last = 0;
*result = 0.;
*abserr = 0.;
if (*limit < 1 || *lenw < *limit << 2) return;
l1 = *limit + 1;
l2 = *limit + l1;
l3 = *limit + l2;
rdqagie(IS, bound, inf, epsabs, epsrel, limit, result, abserr, neval, ier,
&work[1], &work[l1], &work[l2], &work[l3], &iwork[1], last);
return;
} /* Rdqagi */
static
void rdqagie(IntStruct IS, double *bound, int *inf, double *
epsabs, double *epsrel, int *limit, double *result,
double *abserr, int *neval, int *ier, double *alist__,
double *blist, double *rlist, double *elist, int *
iord, int *last)
{
/* System generated locals */
int i__1, i__2;
double d__1, d__2;
/* Local variables */
double area, dres;
int ksgn;
double boun;
int nres;
double area1, area2, area12;
int k;
double small = 0.0, erro12;
int ierro;
double a1, a2, b1, b2, defab1, defab2, oflow;
int ktmin, nrmax;
double uflow;
Rboolean noext;
int iroff1, iroff2, iroff3;
double res3la[3], error1, error2;
int id;
double rlist2[52];
int numrl2;
double defabs, epmach, erlarg = 0.0, abseps, correc = 0.0, errbnd, resabs;
int jupbnd;
double erlast, errmax;
int maxerr;
double reseps;
Rboolean extrap;
double ertest = 0.0, errsum;
/**begin prologue dqagie
***date written 800101 (yymmdd)
***revision date 830518 (yymmdd)
***category no. h2a3a1,h2a4a1
***keywords automatic integrator, infinite intervals,
general-purpose, transformation, extrapolation,
globally adaptive
***author piessens,robert,appl. math & progr. div - k.u.leuven
de doncker,elise,appl. math & progr. div - k.u.leuven
***purpose the routine calculates an approximation result to a given
integral i = integral of f over (bound,+infinity)
or i = integral of f over (-infinity,bound)
or i = integral of f over (-infinity,+infinity),
hopefully satisfying following claim for accuracy
abs(i-result).le.max(epsabs,epsrel*abs(i))
***description
integration over infinite intervals
standard fortran subroutine
f - double precision
function subprogram defining the integrand
function f(x). the actual name for f needs to be
declared e x t e r n a l in the driver program.
bound - double precision
finite bound of integration range
(has no meaning if interval is doubly-infinite)
inf - double precision
indicating the kind of integration range involved
inf = 1 corresponds to (bound,+infinity),
inf = -1 to (-infinity,bound),
inf = 2 to (-infinity,+infinity).
epsabs - double precision
absolute accuracy requested
epsrel - double precision
relative accuracy requested
if epsabs.le.0
and epsrel.lt.max(50*rel.mach.acc.,0.5d-28),
the routine will end with ier = 6.
limit - int
gives an upper bound on the number of subintervals
in the partition of (a,b), limit.ge.1
on return
result - double precision
approximation to the integral
abserr - double precision
estimate of the modulus of the absolute error,
which should equal or exceed abs(i-result)
neval - int
number of integrand evaluations
ier - int
ier = 0 normal and reliable termination of the
routine. it is assumed that the requested
accuracy has been achieved.
- ier.gt.0 abnormal termination of the routine. the
estimates for result and error are less
reliable. it is assumed that the requested
accuracy has not been achieved.
error messages
ier = 1 maximum number of subdivisions allowed
has been achieved. one can allow more
subdivisions by increasing the value of
limit (and taking the according dimension
adjustments into account). however,if
this yields no improvement it is advised
to analyze the integrand in order to
determine the integration difficulties.
if the position of a local difficulty can
be determined (e.g. singularity,
discontinuity within the interval) one
will probably gain from splitting up the
interval at this point and calling the
integrator on the subranges. if possible,
an appropriate special-purpose integrator
should be used, which is designed for
handling the type of difficulty involved.
= 2 the occurrence of roundoff error is
detected, which prevents the requested
tolerance from being achieved.
the error may be under-estimated.
= 3 extremely bad integrand behaviour occurs
at some points of the integration
interval.
= 4 the algorithm does not converge.
roundoff error is detected in the
extrapolation table.
it is assumed that the requested tolerance
cannot be achieved, and that the returned
result is the best which can be obtained.
= 5 the integral is probably divergent, or
slowly convergent. it must be noted that
divergence can occur with any other value
of ier.
= 6 the input is invalid, because
(epsabs.le.0 and
epsrel.lt.max(50*rel.mach.acc.,0.5d-28),
result, abserr, neval, last, rlist(1),
elist(1) and iord(1) are set to zero.
alist(1) and blist(1) are set to 0
and 1 respectively.
alist - double precision
vector of dimension at least limit, the first
last elements of which are the left
end points of the subintervals in the partition
of the transformed integration range (0,1).
blist - double precision
vector of dimension at least limit, the first
last elements of which are the right
end points of the subintervals in the partition
of the transformed integration range (0,1).
rlist - double precision
vector of dimension at least limit, the first
last elements of which are the integral
approximations on the subintervals
elist - double precision
vector of dimension at least limit, the first
last elements of which are the moduli of the
absolute error estimates on the subintervals
iord - int
vector of dimension limit, the first k
elements of which are pointers to the
error estimates over the subintervals,
such that elist(iord(1)), ..., elist(iord(k))
form a decreasing sequence, with k = last
if last.le.(limit/2+2), and k = limit+1-last
otherwise
last - int
number of subintervals actually produced
in the subdivision process
***references (none)
***routines called dqelg,dqk15i,dqpsrt
***end prologue dqagie
the dimension of rlist2 is determined by the value of
limexp in subroutine dqelg.
list of major variables
-----------------------
alist - list of left end points of all subintervals
considered up to now
blist - list of right end points of all subintervals
considered up to now
rlist(i) - approximation to the integral over
(alist(i),blist(i))
rlist2 - array of dimension at least (limexp+2),
containing the part of the epsilon table
wich is still needed for further computations
elist(i) - error estimate applying to rlist(i)
maxerr - pointer to the interval with largest error
estimate
errmax - elist(maxerr)
erlast - error on the interval currently subdivided
(before that subdivision has taken place)
area - sum of the integrals over the subintervals
errsum - sum of the errors over the subintervals
errbnd - requested accuracy max(epsabs,epsrel*
abs(result))
*****1 - variable for the left subinterval
*****2 - variable for the right subinterval
last - index for subdivision
nres - number of calls to the extrapolation routine
numrl2 - number of elements currently in rlist2. if an
appropriate approximation to the compounded
integral has been obtained, it is put in
rlist2(numrl2) after numrl2 has been increased
by one.
small - length of the smallest interval considered up
to now, multiplied by 1.5
erlarg - sum of the errors over the intervals larger
than the smallest interval considered up to now
extrap - logical variable denoting that the routine
is attempting to perform extrapolation. i.e.
before subdividing the smallest interval we
try to decrease the value of erlarg.
noext - logical variable denoting that extrapolation
is no longer allowed (true-value)
machine dependent constants
---------------------------
epmach is the largest relative spacing.
uflow is the smallest positive magnitude.
oflow is the largest positive magnitude. */
/* ***first executable statement dqagie */
/* Parameter adjustments */
--iord;
--elist;
--rlist;
--blist;
--alist__;
/* Function Body */
epmach = DBL_EPSILON;
/* test on validity of parameters */
/* ----------------------------- */
*ier = 0;
*neval = 0;
*last = 0;
*result = 0.;
*abserr = 0.;
alist__[1] = 0.;
blist[1] = 1.;
rlist[1] = 0.;
elist[1] = 0.;
iord[1] = 0;
/* Computing MAX */
d__1 = epmach * 50.;
if (*epsabs <= 0. && (*epsrel < fmax2(epmach * 50., 5e-29))) *ier = 6;
if (*ier == 6) return;
/* first approximation to the integral */
/* ----------------------------------- */
/* determine the interval to be mapped onto (0,1).
if inf = 2 the integral is computed as i = i1+i2, where
i1 = integral of f over (-infinity,0),
i2 = integral of f over (0,+infinity). */
boun = *bound;
if (*inf == 2) {
boun = 0.;
}
rdqk15i(IS, &boun, inf, &c_b6, &c_b7, result, abserr, &defabs, &resabs);
/* test on accuracy */
*last = 1;
rlist[1] = *result;
elist[1] = *abserr;
iord[1] = 1;
dres = fabs(*result);
/* Computing MAX */
d__1 = *epsabs, d__2 = *epsrel * dres;
errbnd = fmax2(d__1,d__2);
if (*abserr <= epmach * 100. * defabs && *abserr > errbnd) *ier = 2;
if (*limit == 1) *ier = 1;
if (*ier != 0 || (*abserr <= errbnd && *abserr != resabs)
|| *abserr == 0.) goto L130;
/* initialization */
/* -------------- */
uflow = DBL_MIN;
oflow = DBL_MAX;
rlist2[0] = *result;
errmax = *abserr;
maxerr = 1;
area = *result;
errsum = *abserr;
*abserr = oflow;
nrmax = 1;
nres = 0;
ktmin = 0;
numrl2 = 2;
extrap = FALSE;
noext = FALSE;
ierro = 0;
iroff1 = 0;
iroff2 = 0;
iroff3 = 0;
ksgn = -1;
if (dres >= (1. - epmach * 50.) * defabs) {
ksgn = 1;
}
/* main do-loop */
/* ------------ */
i__1 = *limit;
for (*last = 2; *last <= i__1; ++(*last)) {
/* bisect the subinterval with nrmax-th largest error estimate. */
a1 = alist__[maxerr];
b1 = (alist__[maxerr] + blist[maxerr]) * .5;
a2 = b1;
b2 = blist[maxerr];
erlast = errmax;
rdqk15i(IS, &boun, inf, &a1, &b1, &area1, &error1, &resabs, &defab1);
rdqk15i(IS, &boun, inf, &a2, &b2, &area2, &error2, &resabs, &defab2);
/* improve previous approximations to integral */
/* and error and test for accuracy. */
area12 = area1 + area2;
erro12 = error1 + error2;
errsum = errsum + erro12 - errmax;
area = area + area12 - rlist[maxerr];
if (defab1 == error1 || defab2 == error2) {
goto L15;
}
if ((d__1 = rlist[maxerr] - area12, fabs(d__1)) > fabs(area12) * 1e-5 ||
erro12 < errmax * .99) {
goto L10;
}
if (extrap) {
++iroff2;
}
if (! extrap) {
++iroff1;
}
L10:
if (*last > 10 && erro12 > errmax) {
++iroff3;
}
L15:
rlist[maxerr] = area1;
rlist[*last] = area2;
/* Computing MAX */
d__1 = *epsabs, d__2 = *epsrel * fabs(area);
errbnd = fmax2(d__1,d__2);
/* test for roundoff error and eventually set error flag. */
if (iroff1 + iroff2 >= 10 || iroff3 >= 20) {
*ier = 2;
}
if (iroff2 >= 5) {
ierro = 3;
}
/* set error flag in the case that the number of */
/* subintervals equals limit. */
if (*last == *limit) {
*ier = 1;
}
/* set error flag in the case of bad integrand behaviour */
/* at some points of the integration range. */
/* Computing MAX */
d__1 = fabs(a1), d__2 = fabs(b2);
if (fmax2(d__1,d__2) <= (epmach * 100. + 1.) * (fabs(a2) + uflow * 1e3))
{
*ier = 4;
}
/* append the newly-created intervals to the list. */
if (error2 > error1) {
goto L20;
}
alist__[*last] = a2;
blist[maxerr] = b1;
blist[*last] = b2;
elist[maxerr] = error1;
elist[*last] = error2;
goto L30;
L20:
alist__[maxerr] = a2;
alist__[*last] = a1;
blist[*last] = b1;
rlist[maxerr] = area2;
rlist[*last] = area1;
elist[maxerr] = error2;
elist[*last] = error1;
/* call subroutine dqpsrt to maintain the descending ordering */
/* in the list of error estimates and select the subinterval */
/* with nrmax-th largest error estimate (to be bisected next). */
L30:
rdqpsrt(limit, last, &maxerr, &errmax, &elist[1], &iord[1], &nrmax);
if (errsum <= errbnd) {
goto L115;
}
if (*ier != 0) {
goto L100;
}
if (*last == 2) {
goto L80;
}
if (noext) {
goto L90;
}
erlarg -= erlast;
if ((d__1 = b1 - a1, fabs(d__1)) > small) {
erlarg += erro12;
}
if (extrap) {
goto L40;
}
/* test whether the interval to be bisected next is the */
/* smallest interval. */
if ((d__1 = blist[maxerr] - alist__[maxerr], fabs(d__1)) > small) {
goto L90;
}
extrap = TRUE;
nrmax = 2;
L40:
if (ierro == 3 || erlarg <= ertest) {
goto L60;
}
/* the smallest interval has the largest error. */
/* before bisecting decrease the sum of the errors over the */
/* larger intervals (erlarg) and perform extrapolation. */
id = nrmax;
jupbnd = *last;
if (*last > *limit / 2 + 2) {
jupbnd = *limit + 3 - *last;
}
i__2 = jupbnd;
for (k = id; k <= i__2; ++k) {
maxerr = iord[nrmax];
errmax = elist[maxerr];
if ((d__1 = blist[maxerr] - alist__[maxerr], fabs(d__1)) > small) {
goto L90;
}
++nrmax;
/* L50: */
}
/* perform extrapolation. */
L60:
++numrl2;
rlist2[numrl2 - 1] = area;
rdqelg(&numrl2, rlist2, &reseps, &abseps, res3la, &nres);
++ktmin;
if (ktmin > 5 && *abserr < errsum * .001) {
*ier = 5;
}
if (abseps >= *abserr) {
goto L70;
}
ktmin = 0;
*abserr = abseps;
*result = reseps;
correc = erlarg;
/* Computing MAX */
d__1 = *epsabs, d__2 = *epsrel * fabs(reseps);
ertest = fmax2(d__1,d__2);
if (*abserr <= ertest) {
goto L100;
}
/* prepare bisection of the smallest interval. */
L70:
if (numrl2 == 1) {
noext = TRUE;
}
if (*ier == 5) {
goto L100;
}
maxerr = iord[1];
errmax = elist[maxerr];
nrmax = 1;
extrap = FALSE;
small *= .5;
erlarg = errsum;
goto L90;
L80:
small = .375;
erlarg = errsum;
ertest = errbnd;
rlist2[1] = area;
L90:
;
}
/* set final result and error estimate. */
/* ------------------------------------ */
L100:
if (*abserr == oflow) {
goto L115;
}
if (*ier + ierro == 0) {
goto L110;
}
if (ierro == 3) {
*abserr += correc;
}
if (*ier == 0) {
*ier = 3;
}
if (*result != 0. && area != 0.) {
goto L105;
}
if (*abserr > errsum) {
goto L115;
}
if (area == 0.) {
goto L130;
}
goto L110;
L105:
if (*abserr / fabs(*result) > errsum / fabs(area)) {
goto L115;
}
/* test on divergence */
L110:
/* Computing MAX */
d__1 = fabs(*result), d__2 = fabs(area);
if (ksgn == -1 && fmax2(d__1,d__2) <= defabs * .01) {
goto L130;
}
if (.01 > *result / area || *result / area > 100. || errsum > fabs(area)) {
*ier = 6;
}
goto L130;
/* compute global integral sum. */
L115:
*result = 0.;
i__1 = *last;
for (k = 1; k <= i__1; ++k) {
*result += rlist[k];
/* L120: */
}
*abserr = errsum;
L130:
*neval = *last * 30 - 15;
if (*inf == 2) {
*neval <<= 1;
}
if (*ier > 2) {
--(*ier);
}
return;
} /* rdqagie_ */
static
void Rdqags(IntStruct IS, double *a, double *b, double *epsabs,
double *epsrel, double *result, double *abserr, int *
neval, int *ier, int *limit, int *lenw, int *last,
int *iwork, double *work)
{
int l1, l2, l3;
/*
***begin prologue dqags
***date written 800101 (yymmdd)
***revision date 830518 (yymmdd)
***category no. h2a1a1
***keywords automatic integrator, general-purpose,
(end-point) singularities, extrapolation,
globally adaptive
***author piessens,robert,appl. math. & progr. div. - k.u.leuven
de doncker,elise,appl. math. & prog. div. - k.u.leuven
***purpose the routine calculates an approximation result to a given
definite integral i = integral of f over (a,b),
hopefully satisfying following claim for accuracy
abs(i-result).le.max(epsabs,epsrel*abs(i)).
***description
computation of a definite integral
standard fortran subroutine
double precision version
parameters
on entry
f - double precision
function subprogram defining the integrand
function f(x). the actual name for f needs to be
declared e x t e r n a l in the driver program.
a - double precision
lower limit of integration
b - double precision
upper limit of integration
epsabs - double precision
absolute accuracy requested
epsrel - double precision
relative accuracy requested
if epsabs.le.0