summaryrefslogtreecommitdiffstats
path: root/tcllib/modules/math/calculus.tcl
blob: 5667a6c5b1e64934817db3ad04d9eedd89a5451f (plain)
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
1001
1002
1003
1004
1005
1006
1007
1008
1009
1010
1011
1012
1013
1014
1015
1016
1017
1018
1019
1020
1021
1022
1023
1024
1025
1026
1027
1028
1029
1030
1031
1032
1033
1034
1035
1036
1037
1038
1039
1040
1041
1042
1043
1044
1045
1046
1047
1048
1049
1050
1051
1052
1053
1054
1055
1056
1057
1058
1059
1060
1061
1062
1063
1064
1065
1066
1067
1068
1069
1070
1071
1072
1073
1074
1075
1076
1077
1078
1079
1080
1081
1082
1083
1084
1085
1086
1087
1088
1089
1090
1091
1092
1093
1094
1095
1096
1097
1098
1099
1100
1101
1102
1103
1104
1105
1106
1107
1108
1109
1110
1111
1112
1113
1114
1115
1116
1117
1118
1119
1120
1121
1122
1123
1124
1125
1126
1127
1128
1129
1130
1131
1132
1133
1134
1135
1136
1137
1138
1139
1140
1141
1142
1143
1144
1145
1146
1147
1148
1149
1150
1151
1152
1153
1154
1155
1156
1157
1158
1159
1160
1161
1162
1163
1164
1165
1166
1167
1168
1169
1170
1171
1172
1173
1174
1175
1176
1177
1178
1179
1180
1181
1182
1183
1184
1185
1186
1187
1188
1189
1190
1191
1192
1193
1194
1195
1196
1197
1198
1199
1200
1201
1202
1203
1204
1205
1206
1207
1208
1209
1210
1211
1212
1213
1214
1215
1216
1217
1218
1219
1220
1221
1222
1223
1224
1225
1226
1227
1228
1229
1230
1231
1232
1233
1234
1235
1236
1237
1238
1239
1240
1241
1242
1243
1244
1245
1246
1247
1248
1249
1250
1251
1252
1253
1254
1255
1256
1257
1258
1259
1260
1261
1262
1263
1264
1265
1266
1267
1268
1269
1270
1271
1272
1273
1274
1275
1276
1277
1278
1279
1280
1281
1282
1283
1284
1285
1286
1287
1288
1289
1290
1291
1292
1293
1294
1295
1296
1297
1298
1299
1300
1301
1302
1303
1304
1305
1306
1307
1308
1309
1310
1311
1312
1313
1314
1315
1316
1317
1318
1319
1320
1321
1322
1323
1324
1325
1326
1327
1328
1329
1330
1331
1332
1333
1334
1335
1336
1337
1338
1339
1340
1341
1342
1343
1344
1345
1346
1347
1348
1349
1350
1351
1352
1353
1354
1355
1356
1357
1358
1359
1360
1361
1362
1363
1364
1365
1366
1367
1368
1369
1370
1371
1372
1373
1374
1375
1376
1377
1378
1379
1380
1381
1382
1383
1384
1385
1386
1387
1388
1389
1390
1391
1392
1393
1394
1395
1396
1397
1398
1399
1400
1401
1402
1403
1404
1405
1406
1407
1408
1409
1410
1411
1412
1413
1414
1415
1416
1417
1418
1419
1420
1421
1422
1423
1424
1425
1426
1427
1428
1429
1430
1431
1432
1433
1434
1435
1436
1437
1438
1439
1440
1441
1442
1443
1444
1445
1446
1447
1448
1449
1450
1451
1452
1453
1454
1455
1456
1457
1458
1459
1460
1461
1462
1463
1464
1465
1466
1467
1468
1469
1470
1471
1472
1473
1474
1475
1476
1477
1478
1479
1480
1481
1482
1483
1484
1485
1486
1487
1488
1489
1490
1491
1492
1493
1494
1495
1496
1497
1498
1499
1500
1501
1502
1503
1504
1505
1506
1507
1508
1509
1510
1511
1512
1513
1514
1515
1516
1517
1518
1519
1520
1521
1522
1523
1524
1525
1526
1527
1528
1529
1530
1531
1532
1533
1534
1535
1536
1537
1538
1539
1540
1541
1542
1543
1544
1545
1546
1547
1548
1549
1550
1551
1552
1553
1554
1555
1556
1557
1558
1559
1560
1561
1562
1563
1564
1565
1566
1567
1568
1569
1570
1571
1572
1573
1574
1575
1576
1577
1578
1579
1580
1581
1582
1583
1584
1585
1586
1587
1588
1589
1590
1591
1592
1593
1594
1595
1596
1597
1598
1599
1600
1601
1602
1603
1604
1605
1606
1607
1608
1609
1610
1611
1612
1613
1614
1615
1616
1617
1618
1619
1620
1621
1622
1623
1624
1625
1626
1627
1628
1629
1630
1631
1632
1633
1634
1635
1636
1637
1638
1639
1640
1641
1642
1643
1644
1645
# calculus.tcl --
#    Package that implements several basic numerical methods, such
#    as the integration of a one-dimensional function and the
#    solution of a system of first-order differential equations.
#
# Copyright (c) 2002, 2003, 2004, 2006 by Arjen Markus.
# Copyright (c) 2004 by Kevin B. Kenny.  All rights reserved.
# See the file "license.terms" for information on usage and redistribution
# of this file, and for a DISCLAIMER OF ALL WARRANTIES.
#
# RCS: @(#) $Id: calculus.tcl,v 1.15 2008/10/08 03:30:48 andreas_kupries Exp $

package require Tcl 8.4
package require math::interpolate
package provide math::calculus 0.8.1

# math::calculus --
#    Namespace for the commands

namespace eval ::math::calculus {

    namespace import ::math::interpolate::neville

    namespace import ::math::expectDouble ::math::expectInteger

    namespace export \
	integral integralExpr integral2D integral3D \
	qk15 qk15_detailed \
	eulerStep heunStep rungeKuttaStep           \
	boundaryValueSecondOrder solveTriDiagonal   \
	newtonRaphson newtonRaphsonParameters
    namespace export \
        integral2D_2accurate integral3D_accurate

    namespace export romberg romberg_infinity
    namespace export romberg_sqrtSingLower romberg_sqrtSingUpper
    namespace export romberg_powerLawLower romberg_powerLawUpper
    namespace export romberg_expLower romberg_expUpper

    namespace export regula_falsi

    variable nr_maxiter    20
    variable nr_tolerance   0.001

}

# integral --
#    Integrate a function over a given interval using the Simpson rule
#
# Arguments:
#    begin       Start of the interval
#    end         End of the interval
#    nosteps     Number of steps in which to divide the interval
#    func        Name of the function to be integrated (takes one
#                argument)
# Return value:
#    Computed integral
#
proc ::math::calculus::integral { begin end nosteps func } {

   set delta    [expr {($end-$begin)/double($nosteps)}]
   set hdelta   [expr {$delta/2.0}]
   set result   0.0
   set xval     $begin
   set func_end [uplevel 1 $func $xval]
   for { set i 1 } { $i <= $nosteps } { incr i } {
      set func_begin  $func_end
      set func_middle [uplevel 1 $func [expr {$xval+$hdelta}]]
      set func_end    [uplevel 1 $func [expr {$xval+$delta}]]
      set result      [expr  {$result+$func_begin+4.0*$func_middle+$func_end}]

      set xval        [expr {$begin+double($i)*$delta}]
   }

   return [expr {$result*$delta/6.0}]
}

# integralExpr --
#    Integrate an expression with "x" as the integrate according to the
#    Simpson rule
#
# Arguments:
#    begin       Start of the interval
#    end         End of the interval
#    nosteps     Number of steps in which to divide the interval
#    expression  Expression with "x" as the integrate
# Return value:
#    Computed integral
#
proc ::math::calculus::integralExpr { begin end nosteps expression } {

   set delta    [expr {($end-$begin)/double($nosteps)}]
   set hdelta   [expr {$delta/2.0}]
   set result   0.0
   set x        $begin
   # FRINK: nocheck
   set func_end [expr $expression]
   for { set i 1 } { $i <= $nosteps } { incr i } {
      set func_begin  $func_end
      set x           [expr {$x+$hdelta}]
       # FRINK: nocheck
      set func_middle [expr $expression]
      set x           [expr {$x+$hdelta}]
       # FRINK: nocheck
      set func_end    [expr $expression]
      set result      [expr {$result+$func_begin+4.0*$func_middle+$func_end}]

      set x           [expr {$begin+double($i)*$delta}]
   }

   return [expr {$result*$delta/6.0}]
}

# integral2D --
#    Integrate a given fucntion of two variables over a block,
#    using bilinear interpolation (for this moment: block function)
#
# Arguments:
#    xinterval   Start, stop and number of steps of the "x" interval
#    yinterval   Start, stop and number of steps of the "y" interval
#    func        Function of the two variables to be integrated
# Return value:
#    Computed integral
#
proc ::math::calculus::integral2D { xinterval yinterval func } {

   foreach { xbegin xend xnumber } $xinterval { break }
   foreach { ybegin yend ynumber } $yinterval { break }

   set xdelta    [expr {($xend-$xbegin)/double($xnumber)}]
   set ydelta    [expr {($yend-$ybegin)/double($ynumber)}]
   set hxdelta   [expr {$xdelta/2.0}]
   set hydelta   [expr {$ydelta/2.0}]
   set result   0.0
   set dxdy      [expr {$xdelta*$ydelta}]
   for { set j 0 } { $j < $ynumber } { incr j } {
      set y [expr {$ybegin+$hydelta+double($j)*$ydelta}]
      for { set i 0 } { $i < $xnumber } { incr i } {
         set x           [expr {$xbegin+$hxdelta+double($i)*$xdelta}]
         set func_value  [uplevel 1 $func $x $y]
         set result      [expr {$result+$func_value}]
      }
   }

   return [expr {$result*$dxdy}]
}

# integral3D --
#    Integrate a given fucntion of two variables over a block,
#    using trilinear interpolation (for this moment: block function)
#
# Arguments:
#    xinterval   Start, stop and number of steps of the "x" interval
#    yinterval   Start, stop and number of steps of the "y" interval
#    zinterval   Start, stop and number of steps of the "z" interval
#    func        Function of the three variables to be integrated
# Return value:
#    Computed integral
#
proc ::math::calculus::integral3D { xinterval yinterval zinterval func } {

   foreach { xbegin xend xnumber } $xinterval { break }
   foreach { ybegin yend ynumber } $yinterval { break }
   foreach { zbegin zend znumber } $zinterval { break }

   set xdelta    [expr {($xend-$xbegin)/double($xnumber)}]
   set ydelta    [expr {($yend-$ybegin)/double($ynumber)}]
   set zdelta    [expr {($zend-$zbegin)/double($znumber)}]
   set hxdelta   [expr {$xdelta/2.0}]
   set hydelta   [expr {$ydelta/2.0}]
   set hzdelta   [expr {$zdelta/2.0}]
   set result   0.0
   set dxdydz    [expr {$xdelta*$ydelta*$zdelta}]
   for { set k 0 } { $k < $znumber } { incr k } {
      set z [expr {$zbegin+$hzdelta+double($k)*$zdelta}]
      for { set j 0 } { $j < $ynumber } { incr j } {
         set y [expr {$ybegin+$hydelta+double($j)*$ydelta}]
         for { set i 0 } { $i < $xnumber } { incr i } {
            set x           [expr {$xbegin+$hxdelta+double($i)*$xdelta}]
            set func_value  [uplevel 1 $func $x $y $z]
            set result      [expr {$result+$func_value}]
         }
      }
   }

   return [expr {$result*$dxdydz}]
}

# integral2D_accurate --
#    Integrate a given function of two variables over a block,
#    using a four-point quadrature formula
#
# Arguments:
#    xinterval   Start, stop and number of steps of the "x" interval
#    yinterval   Start, stop and number of steps of the "y" interval
#    func        Function of the two variables to be integrated
# Return value:
#    Computed integral
#
proc ::math::calculus::integral2D_accurate { xinterval yinterval func } {

    foreach { xbegin xend xnumber } $xinterval { break }
    foreach { ybegin yend ynumber } $yinterval { break }

    set alpha     [expr {sqrt(2.0/3.0)}]
    set minalpha  [expr {-$alpha}]
    set dpoints   [list $alpha 0.0 $minalpha 0.0 0.0 $alpha 0.0 $minalpha]

    set xdelta    [expr {($xend-$xbegin)/double($xnumber)}]
    set ydelta    [expr {($yend-$ybegin)/double($ynumber)}]
    set hxdelta   [expr {$xdelta/2.0}]
    set hydelta   [expr {$ydelta/2.0}]
    set result   0.0
    set dxdy      [expr {0.25*$xdelta*$ydelta}]

    for { set j 0 } { $j < $ynumber } { incr j } {
        set y [expr {$ybegin+$hydelta+double($j)*$ydelta}]
        for { set i 0 } { $i < $xnumber } { incr i } {
            set x [expr {$xbegin+$hxdelta+double($i)*$xdelta}]

            foreach {dx dy} $dpoints {
                set x1          [expr {$x+$dx}]
                set y1          [expr {$y+$dy}]
                set func_value  [uplevel 1 $func $x1 $y1]
                set result      [expr {$result+$func_value}]
            }
        }
    }

    return [expr {$result*$dxdy}]
}

# integral3D_accurate --
#    Integrate a given function of three variables over a block,
#    using an 8-point quadrature formula
#
# Arguments:
#    xinterval   Start, stop and number of steps of the "x" interval
#    yinterval   Start, stop and number of steps of the "y" interval
#    zinterval   Start, stop and number of steps of the "z" interval
#    func        Function of the three variables to be integrated
# Return value:
#    Computed integral
#
proc ::math::calculus::integral3D_accurate { xinterval yinterval zinterval func } {

    foreach { xbegin xend xnumber } $xinterval { break }
    foreach { ybegin yend ynumber } $yinterval { break }
    foreach { zbegin zend znumber } $zinterval { break }

    set alpha     [expr {sqrt(1.0/3.0)}]
    set minalpha  [expr {-$alpha}]

    set dpoints   [list $alpha    $alpha    $alpha    \
                        $alpha    $alpha    $minalpha \
                        $alpha    $minalpha $alpha    \
                        $alpha    $minalpha $minalpha \
                        $minalpha $alpha    $alpha    \
                        $minalpha $alpha    $minalpha \
                        $minalpha $minalpha $alpha    \
                        $minalpha $minalpha $minalpha ]

    set xdelta    [expr {($xend-$xbegin)/double($xnumber)}]
    set ydelta    [expr {($yend-$ybegin)/double($ynumber)}]
    set zdelta    [expr {($zend-$zbegin)/double($znumber)}]
    set hxdelta   [expr {$xdelta/2.0}]
    set hydelta   [expr {$ydelta/2.0}]
    set hzdelta   [expr {$zdelta/2.0}]
    set result    0.0
    set dxdydz    [expr {0.125*$xdelta*$ydelta*$zdelta}]

    for { set k 0 } { $k < $znumber } { incr k } {
        set z [expr {$zbegin+$hzdelta+double($k)*$zdelta}]
        for { set j 0 } { $j < $ynumber } { incr j } {
            set y [expr {$ybegin+$hydelta+double($j)*$ydelta}]
            for { set i 0 } { $i < $xnumber } { incr i } {
                set x [expr {$xbegin+$hxdelta+double($i)*$xdelta}]

                foreach {dx dy dz} $dpoints {
                    set x1 [expr {$x+$dx}]
                    set y1 [expr {$y+$dy}]
                    set z1 [expr {$z+$dz}]
                    set func_value  [uplevel 1 $func $x1 $y1 $z1]
                    set result      [expr {$result+$func_value}]
                }
            }
        }
    }

    return [expr {$result*$dxdydz}]
}

# eulerStep --
#    Integrate a system of ordinary differential equations of the type
#    x' = f(x,t), where x is a vector of quantities. Integration is
#    done over a single step according to Euler's method.
#
# Arguments:
#    t           Start value of independent variable (time for instance)
#    tstep       Step size of interval
#    xvec        Vector of dependent values at the start
#    func        Function taking the arguments t and xvec to return
#                the derivative of each dependent variable.
# Return value:
#    List of values at the end of the step
#
proc ::math::calculus::eulerStep { t tstep xvec func } {

   set xderiv   [uplevel 1 $func $t [list $xvec]]
   set result   {}
   foreach xv $xvec dx $xderiv {
      set xnew [expr {$xv+$tstep*$dx}]
      lappend result $xnew
   }

   return $result
}

# heunStep --
#    Integrate a system of ordinary differential equations of the type
#    x' = f(x,t), where x is a vector of quantities. Integration is
#    done over a single step according to Heun's method.
#
# Arguments:
#    t           Start value of independent variable (time for instance)
#    tstep       Step size of interval
#    xvec        Vector of dependent values at the start
#    func        Function taking the arguments t and xvec to return
#                the derivative of each dependent variable.
# Return value:
#    List of values at the end of the step
#
proc ::math::calculus::heunStep { t tstep xvec func } {

   #
   # Predictor step
   #
   set funcq    [uplevel 1 namespace which -command $func]
   set xpred    [eulerStep $t $tstep $xvec $funcq]

   #
   # Corrector step
   #
   set tcorr    [expr {$t+$tstep}]
   set xcorr    [eulerStep $tcorr $tstep $xpred $funcq]

   set result   {}
   foreach xv $xvec xc $xcorr {
      set xnew [expr {0.5*($xv+$xc)}]
      lappend result $xnew
   }

   return $result
}

# rungeKuttaStep --
#    Integrate a system of ordinary differential equations of the type
#    x' = f(x,t), where x is a vector of quantities. Integration is
#    done over a single step according to Runge-Kutta 4th order.
#
# Arguments:
#    t           Start value of independent variable (time for instance)
#    tstep       Step size of interval
#    xvec        Vector of dependent values at the start
#    func        Function taking the arguments t and xvec to return
#                the derivative of each dependent variable.
# Return value:
#    List of values at the end of the step
#
proc ::math::calculus::rungeKuttaStep { t tstep xvec func } {

   set funcq    [uplevel 1 namespace which -command $func]

   #
   # Four steps:
   # - k1 = tstep*func(t,x0)
   # - k2 = tstep*func(t+0.5*tstep,x0+0.5*k1)
   # - k3 = tstep*func(t+0.5*tstep,x0+0.5*k2)
   # - k4 = tstep*func(t+    tstep,x0+    k3)
   # - x1 = x0 + (k1+2*k2+2*k3+k4)/6
   #
   set tstep2   [expr {$tstep/2.0}]
   set tstep6   [expr {$tstep/6.0}]

   set xk1      [$funcq $t $xvec]
   set xvec2    {}
   foreach x1 $xvec xv $xk1 {
      lappend xvec2 [expr {$x1+$tstep2*$xv}]
   }
   set xk2      [$funcq [expr {$t+$tstep2}] $xvec2]

   set xvec3    {}
   foreach x1 $xvec xv $xk2 {
      lappend xvec3 [expr {$x1+$tstep2*$xv}]
   }
   set xk3      [$funcq [expr {$t+$tstep2}] $xvec3]

   set xvec4    {}
   foreach x1 $xvec xv $xk3 {
      lappend xvec4 [expr {$x1+$tstep*$xv}]
   }
   set xk4      [$funcq [expr {$t+$tstep}] $xvec4]

   set result   {}
   foreach x0 $xvec k1 $xk1 k2 $xk2 k3 $xk3 k4 $xk4 {
      set dx [expr {$k1+2.0*$k2+2.0*$k3+$k4}]
      lappend result [expr {$x0+$dx*$tstep6}]
   }

   return $result
}

# boundaryValueSecondOrder --
#    Integrate a second-order differential equation and solve for
#    given boundary values.
#
#    The equation is (see the documentation):
#       d       dy   d
#       -- A(x) -- + -- B(x) y + C(x) y = D(x)
#       dx      dx   dx
#
#    The procedure uses finite differences and tridiagonal matrices to
#    solve the equation. The boundary values are put in the matrix
#    directly.
#
# Arguments:
#    coeff_func  Name of triple-valued function for coefficients A, B, C
#    force_func  Name of the function providing the force term D(x)
#    leftbnd     Left boundary condition (list of: xvalue, boundary
#                value or keyword zero-flux, zero-derivative)
#    rightbnd    Right boundary condition (ditto)
#    nostep      Number of steps
# Return value:
#    List of x-values and calculated values (x1, y1, x2, y2, ...)
#
proc ::math::calculus::boundaryValueSecondOrder {
   coeff_func force_func leftbnd rightbnd nostep } {

   set coeffq    [uplevel 1 namespace which -command $coeff_func]
   set forceq    [uplevel 1 namespace which -command $force_func]

   if { [llength $leftbnd] != 2 || [llength $rightbnd] != 2 } {
      error "Boundary condition(s) incorrect"
   }
   if { $nostep < 1 } {
      error "Number of steps must be larger/equal 1"
   }

   #
   # Set up the matrix, as three different lists and the
   # righthand side as the fourth
   #
   set xleft  [lindex $leftbnd 0]
   set xright [lindex $rightbnd 0]
   set xstep  [expr {($xright-$xleft)/double($nostep)}]

   set acoeff {}
   set bcoeff {}
   set ccoeff {}
   set dvalue {}

   set x $xleft
   foreach {A B C} [$coeffq $x] { break }

   set A1 [expr {$A/$xstep-0.5*$B}]
   set B1 [expr {$A/$xstep+0.5*$B+0.5*$C*$xstep}]
   set C1 0.0

   for { set i 1 } { $i <= $nostep } { incr i } {
      set x [expr {$xleft+double($i)*$xstep}]
      if { [expr {abs($x)-0.5*abs($xstep)}] < 0.0 } {
         set x 0.0
      }
      foreach {A B C} [$coeffq $x] { break }

      set A2 0.0
      set B2 [expr {$A/$xstep-0.5*$B+0.5*$C*$xstep}]
      set C2 [expr {$A/$xstep+0.5*$B}]
      lappend acoeff [expr {$A1+$A2}]
      lappend bcoeff [expr {-$B1-$B2}]
      lappend ccoeff [expr {$C1+$C2}]
      set A1 [expr {$A/$xstep-0.5*$B}]
      set B1 [expr {$A/$xstep+0.5*$B+0.5*$C*$xstep}]
      set C1 0.0
   }
   set xvec {}
   for { set i 0 } { $i < $nostep } { incr i } {
      set x [expr {$xleft+(0.5+double($i))*$xstep}]
      if { [expr {abs($x)-0.25*abs($xstep)}] < 0.0 } {
         set x 0.0
      }
      lappend xvec   $x
      lappend dvalue [expr {$xstep*[$forceq $x]}]
   }

   #
   # Substitute the boundary values
   #
   set A  [lindex $acoeff 0]
   set D  [lindex $dvalue 0]
   set D1 [expr {$D-$A*[lindex $leftbnd 1]}]
   set C  [lindex $ccoeff end]
   set D  [lindex $dvalue end]
   set D2 [expr {$D-$C*[lindex $rightbnd 1]}]
   set dvalue [concat $D1 [lrange $dvalue 1 end-1] $D2]

   set yvec [solveTriDiagonal [lrange $acoeff 1 end] $bcoeff [lrange $ccoeff 0 end-1] $dvalue]

   foreach x $xvec y $yvec {
      lappend result $x $y
   }
   return $result
}

# solveTriDiagonal --
#    Solve a system of equations Ax = b where A is a tridiagonal matrix
#
# Arguments:
#    acoeff      Values on lower diagonal
#    bcoeff      Values on main diagonal
#    ccoeff      Values on upper diagonal
#    dvalue      Values on righthand side
# Return value:
#    List of values forming the solution
#
proc ::math::calculus::solveTriDiagonal { acoeff bcoeff ccoeff dvalue } {

   set nostep [llength $acoeff]
   #
   # First step: Gauss-elimination
   #
   set B [lindex $bcoeff 0]
   set C [lindex $ccoeff 0]
   set D [lindex $dvalue 0]
   set acoeff  [concat 0.0 $acoeff]
   set bcoeff2 [list $B]
   set dvalue2 [list $D]
   for { set i 1 } { $i <= $nostep } { incr i } {
      set A2    [lindex $acoeff $i]
      set B2    [lindex $bcoeff $i]
      set D2    [lindex $dvalue $i]
      set ratab [expr {$A2/double($B)}]
      set B2    [expr {$B2-$ratab*$C}]
      set D2    [expr {$D2-$ratab*$D}]
      lappend bcoeff2 $B2
      lappend dvalue2 $D2
      set B     $B2
      set C     [lindex $ccoeff $i]
      set D     $D2
   }

   #
   # Second step: substitution
   #
   set yvec {}
   set B [lindex $bcoeff2 end]
   set D [lindex $dvalue2 end]
   set y [expr {$D/$B}]
   for { set i [expr {$nostep-1}] } { $i >= 0 } { incr i -1 } {
      set yvec  [concat $y $yvec]
      set B     [lindex $bcoeff2 $i]
      set C     [lindex $ccoeff  $i]
      set D     [lindex $dvalue2 $i]
      set y     [expr {($D-$C*$y)/$B}]
   }
   set yvec [concat $y $yvec]

   return $yvec
}

# newtonRaphson --
#    Determine the root of an equation via the Newton-Raphson method
#
# Arguments:
#    func        Function (proc) in x
#    deriv       Derivative (proc) of func w.r.t. x
#    initval     Initial value for x
# Return value:
#    Estimate of root
#
proc ::math::calculus::newtonRaphson { func deriv initval } {
   variable nr_maxiter
   variable nr_tolerance

   set funcq  [uplevel 1 namespace which -command $func]
   set derivq [uplevel 1 namespace which -command $deriv]

   set value $initval
   set diff  [expr {10.0*$nr_tolerance}]

   for { set i 0 } { $i < $nr_maxiter } { incr i } {
      if { $diff < $nr_tolerance } {
         break
      }

      set newval [expr {$value-[$funcq $value]/[$derivq $value]}]
      if { $value != 0.0 } {
         set diff   [expr {abs($newval-$value)/abs($value)}]
      } else {
         set diff   [expr {abs($newval-$value)}]
      }
      set value $newval
   }

   return $newval
}

# newtonRaphsonParameters --
#    Set the parameters for the Newton-Raphson method
#
# Arguments:
#    maxiter     Maximum number of iterations
#    tolerance   Relative precisiion of the result
# Return value:
#    None
#
proc ::math::calculus::newtonRaphsonParameters { maxiter tolerance } {
   variable nr_maxiter
   variable nr_tolerance

   if { $maxiter > 0 } {
      set nr_maxiter $maxiter
   }
   if { $tolerance > 0 } {
      set nr_tolerance $tolerance
   }
}

#----------------------------------------------------------------------
#
# midpoint --
#
#	Perform one set of steps in evaluating an integral using the
#	midpoint method.
#
# Usage:
#	midpoint f a b s ?n?
#
# Parameters:
#	f - function to integrate
#	a - One limit of integration
#	b - Other limit of integration.  a and b need not be in ascending
#	    order.
#	s - Value returned from a previous call to midpoint (see below)
#	n - Step number (see below)
#
# Results:
#	Returns an estimate of the integral obtained by dividing the
#	interval into 3**n equal intervals and using the midpoint rule.
#
# Side effects:
#	f is evaluated 2*3**(n-1) times and may have side effects.
#
# The 'midpoint' procedure is designed for successive approximations.
# It should be called initially with n==0.  On this initial call, s
# is ignored.  The function is evaluated at the midpoint of the interval, and
# the value is multiplied by the width of the interval to give the
# coarsest possible estimate of the integral.
#
# On each iteration except the first, n should be incremented by one,
# and the previous value returned from [midpoint] should be supplied
# as 's'.  The function will be evaluated at additional points
# to give a total of 3**n equally spaced points, and the estimate
# of the integral will be updated and returned
#
# Under normal circumstances, user code will not call this function
# directly. Instead, it will use ::math::calculus::romberg to
# do error control and extrapolation to a zero step size.
#
#----------------------------------------------------------------------

proc ::math::calculus::midpoint { f a b { n 0 } { s 0. } } {

    if { $n == 0 } {

	# First iteration.  Simply evaluate the function at the midpoint
	# of the interval.

	set cmd $f; lappend cmd [expr { 0.5 * ( $a + $b ) }]; set v [eval $cmd]
	return [expr { ( $b - $a ) * $v }]

    } else {

	# Subsequent iterations. We've divided the interval into
	# $it subintervals.  Evaluate the function at the 1/3 and
	# 2/3 points of each subinterval.  Then update the estimate
	# of the integral that we produced on the last step with
	# the new sum.

	set it [expr { pow( 3, $n-1 ) }]
	set h [expr { ( $b - $a ) / ( 3. * $it ) }]
	set h2 [expr { $h + $h }]
	set x [expr { $a + 0.5 * $h }]
	set sum 0
	for { set j 0 } { $j < $it } { incr j } {
	    set cmd $f; lappend cmd $x; set y [eval $cmd]
	    set sum [expr { $sum + $y }]
	    set x [expr { $x + $h2 }]
	    set cmd $f; lappend cmd $x; set y [eval $cmd]
	    set sum [expr { $sum + $y }]
	    set x [expr { $x + $h}]
	}
	return [expr { ( $s + ( $b - $a ) * $sum / $it ) / 3. }]

    }
}

#----------------------------------------------------------------------
#
# romberg --
#	
#	Compute the integral of a function over an interval using
#	Romberg's method.
#
# Usage:
#	romberg f a b ?-option value?...
#
# Parameters:
#	f - Function to integrate.  Must be a single Tcl command,
#	    to which will be appended the abscissa at which the function
#	    should be evaluated.  f should be analytic over the
#	    region of integration, but may have a removable singularity
#	    at either endpoint.
#	a - One bound of the interval
#	b - The other bound of the interval.  a and b need not be in
#	    ascending order.
#
# Options:
#	-abserror ABSERROR
#		Requests that the integration be performed to make
#		the estimated absolute error of the integral less than
#		the given value.  Default is 1.e-10.
#	-relerror RELERROR
#		Requests that the integration be performed to make
#		the estimated absolute error of the integral less than
#		the given value.  Default is 1.e-6.
#	-degree N
#		Specifies the degree of the polynomial that will be
#		used to extrapolate to a zero step size.  -degree 0
#		requests integration with the midpoint rule; -degree 1
#		is equivalent to Simpson's 3/8 rule; higher degrees
#		are difficult to describe but (within reason) give
#		faster convergence for smooth functions.  Default is
#		-degree 4.
#	-maxiter N
#		Specifies the maximum number of triplings of the
#		number of steps to take in integration.  At most
#		3**N function evaluations will be performed in
#		integrating with -maxiter N.  The integration
#		will terminate at that time, even if the result
#		satisfies neither the -relerror nor -abserror tests.
#
# Results:
#	Returns a two-element list.  The first element is the estimated
#	value of the integral; the second is the estimated absolute
#	error of the value.
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg { f a b args } {

    # Replace f with a context-independent version

    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]

    # Assign default parameters

    array set params {
	-abserror 1.0e-10
	-degree 4
	-relerror 1.0e-6
	-maxiter 14
    }

    # Extract parameters

    if { ( [llength $args] % 2 ) != 0 } {
        return -code error -errorcode [list romberg wrongNumArgs] \
            "wrong \# args, should be\
                 \"[lreplace [info level 0] 1 end \
                         f x1 x2 ?-option value?...]\""
    }
    foreach { key value } $args {
        if { ![info exists params($key)] } {
            return -code error -errorcode [list romberg badoption $key] \
                "unknown option \"$key\",\
                     should be -abserror, -degree, -relerror, or -maxiter"
        }
        set params($key) $value
    }

    # Check params

    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { ![string is double -strict $params(-abserror)] } {
	return -code error [expectDouble $params(-abserror)]
    }
    if { ![string is integer -strict $params(-degree)] } {
	return -code error [expectInteger $params(-degree)]
    }
    if { ![string is integer -strict $params(-maxiter)] } {
	return -code error [expectInteger $params(-maxiter)]
    }
    if { ![string is double -strict $params(-relerror)] } {
	return -code error [expectDouble $params(-relerror)]
    }
    foreach key {-abserror -degree -maxiter -relerror} {
	if { $params($key) <= 0 } {
	    return -code error -errorcode [list romberg notPositive $key] \
		"$key must be positive"
	}
    }
    if { $params(-maxiter) <= $params(-degree) } {
	return -code error -errorcode [list romberg tooFewIter] \
	    "-maxiter must be greater than -degree"
    }

    # Create lists of step size and sum with the given number of steps.

    set x [list]
    set y [list]
    set s 0;				# Current best estimate of integral
    set indx end-$params(-degree)
    set pow3 1.;			# Current step size (times b-a)

    # Perform successive integrations, tripling the number of steps each time

    for { set i 0 } { $i < $params(-maxiter) } { incr i } {
	set s [midpoint $f $a $b $i $s]
	lappend x $pow3
	lappend y $s
	set pow3 [expr { $pow3 / 9. }]

	# Once $degree steps have been done, start Richardson extrapolation
	# to a zero step size.

	if { $i >= $params(-degree) } {
	    set x [lrange $x $indx end]
	    set y [lrange $y $indx end]
	    foreach {estimate err} [neville $x $y 0.] break
	    if { $err < $params(-abserror)
		 || $err < $params(-relerror) * abs($estimate) } {
		return [list $estimate $err]
	    }
	}
    }

    # If -maxiter iterations have been done, give up, and return
    # with the current error estimate.

    return [list $estimate $err]
}

#----------------------------------------------------------------------
#
# u_infinity --
#	Change of variable for integrating over a half-infinite
#	interval
#
# Parameters:
#	f - Function being integrated
#	u - 1/x, where x is the abscissa where f is to be evaluated
#
# Results:
#	Returns f(1/u)/(u**2)
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_infinity { f u } {
    set cmd $f
    lappend cmd [expr { 1.0 / $u }]
    set y [eval $cmd]
    return [expr { $y / ( $u * $u ) }]
}

#----------------------------------------------------------------------
#
# romberg_infinity --
#	Evaluate a function on a half-open interval
#
# Usage:
#	Same as 'romberg'
#
# The 'romberg_infinity' procedure performs Romberg integration on
# an interval [a,b] where an infinite a or b may be represented by
# a large number (e.g. 1.e30).  It operates by a change of variable;
# instead of integrating f(x) from a to b, it makes a change
# of variable u = 1/x, and integrates from 1/b to 1/a f(1/u)/u**2 du.
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_infinity { f a b args } {
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a * $b <= 0. } {
        return -code error -errorcode {romberg_infinity cross-axis} \
            "limits of integration have opposite sign"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set f [list u_infinity $f]
    return [eval [linsert $args 0 \
                      romberg $f [expr { 1.0 / $b }] [expr { 1.0 / $a }]]]
}

#----------------------------------------------------------------------
#
# u_sqrtSingLower --
#	Change of variable for integrating over an interval with
#	an inverse square root singularity at the lower bound.
#
# Parameters:
#	f - Function being integrated
#	a - Lower bound
#	u - sqrt(x-a), where x is the abscissa where f is to be evaluated
#
# Results:
#	Returns 2 * u * f( a + u**2 )
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_sqrtSingLower { f a u } {
    set cmd $f
    lappend cmd [expr { $a + $u * $u }]
    set y [eval $cmd]
    return [expr { 2. * $u * $y }]
}

#----------------------------------------------------------------------
#
# u_sqrtSingUpper --
#	Change of variable for integrating over an interval with
#	an inverse square root singularity at the upper bound.
#
# Parameters:
#	f - Function being integrated
#	b - Upper bound
#	u - sqrt(b-x), where x is the abscissa where f is to be evaluated
#
# Results:
#	Returns 2 * u * f( b - u**2 )
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_sqrtSingUpper { f b u } {
    set cmd $f
    lappend cmd [expr { $b - $u * $u }]
    set y [eval $cmd]
    return [expr { 2. * $u * $y }]
}

#----------------------------------------------------------------------
#
# math::calculus::romberg_sqrtSingLower --
#	Integrate a function with an inverse square root singularity
#	at the lower bound
#
# Usage:
#	Same as 'romberg'
#
# The 'romberg_sqrtSingLower' procedure is a wrapper for 'romberg'
# for integrating a function with an inverse square root singularity
# at the lower bound of the interval.  It works by making the change
# of variable u = sqrt( x-a ).
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_sqrtSingLower { f a b args } {
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a >= $b } {
	return -code error "limits of integration out of order"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set f [list u_sqrtSingLower $f $a]
    return [eval [linsert $args 0 \
			     romberg $f 0 [expr { sqrt( $b - $a ) }]]]
}

#----------------------------------------------------------------------
#
# math::calculus::romberg_sqrtSingUpper --
#	Integrate a function with an inverse square root singularity
#	at the upper bound
#
# Usage:
#	Same as 'romberg'
#
# The 'romberg_sqrtSingUpper' procedure is a wrapper for 'romberg'
# for integrating a function with an inverse square root singularity
# at the upper bound of the interval.  It works by making the change
# of variable u = sqrt( b-x ).
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_sqrtSingUpper { f a b args } {
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a >= $b } {
	return -code error "limits of integration out of order"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set f [list u_sqrtSingUpper $f $b]
    return [eval [linsert $args 0 \
		      romberg $f 0. [expr { sqrt( $b - $a ) }]]]
}

#----------------------------------------------------------------------
#
# u_powerLawLower --
#	Change of variable for integrating over an interval with
#	an integrable power law singularity at the lower bound.
#
# Parameters:
#	f - Function being integrated
#	gammaover1mgamma - gamma / (1 - gamma), where gamma is the power
#	oneover1mgamma - 1 / (1 - gamma), where gamma is the power
#	a - Lower limit of integration
#	u - Changed variable u == (x-a)**(1-gamma)
#
# Results:
#	Returns u**(1/1-gamma) * f(a + u**(1/1-gamma) ).
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_powerLawLower { f gammaover1mgamma oneover1mgamma
					 a u } {
    set cmd $f
    lappend cmd [expr { $a + pow( $u, $oneover1mgamma ) }]
    set y [eval $cmd]
    return [expr { $y * pow( $u, $gammaover1mgamma ) }]
}

#----------------------------------------------------------------------
#
# math::calculus::romberg_powerLawLower --
#	Integrate a function with an integrable power law singularity
#	at the lower bound
#
# Usage:
#	romberg_powerLawLower gamma f a b ?-option value...?
#
# Parameters:
#	gamma - Power (0<gamma<1) of the singularity
#	f - Function to integrate.  Must be a single Tcl command,
#	    to which will be appended the abscissa at which the function
#	    should be evaluated.  f is expected to have an integrable
#	    power law singularity at the lower endpoint; that is, the
#	    integrand is expected to diverge as (x-a)**gamma.
#	a - One bound of the interval
#	b - The other bound of the interval.  a and b need not be in
#	    ascending order.
#
# Options:
#	-abserror ABSERROR
#		Requests that the integration be performed to make
#		the estimated absolute error of the integral less than
#		the given value.  Default is 1.e-10.
#	-relerror RELERROR
#		Requests that the integration be performed to make
#		the estimated absolute error of the integral less than
#		the given value.  Default is 1.e-6.
#	-degree N
#		Specifies the degree of the polynomial that will be
#		used to extrapolate to a zero step size.  -degree 0
#		requests integration with the midpoint rule; -degree 1
#		is equivalent to Simpson's 3/8 rule; higher degrees
#		are difficult to describe but (within reason) give
#		faster convergence for smooth functions.  Default is
#		-degree 4.
#	-maxiter N
#		Specifies the maximum number of triplings of the
#		number of steps to take in integration.  At most
#		3**N function evaluations will be performed in
#		integrating with -maxiter N.  The integration
#		will terminate at that time, even if the result
#		satisfies neither the -relerror nor -abserror tests.
#
# Results:
#	Returns a two-element list.  The first element is the estimated
#	value of the integral; the second is the estimated absolute
#	error of the value.
#
# The 'romberg_sqrtSingLower' procedure is a wrapper for 'romberg'
# for integrating a function with an integrable power law singularity
# at the lower bound of the interval.  It works by making the change
# of variable u = (x-a)**(1-gamma).
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_powerLawLower { gamma f a b args } {
    if { ![string is double -strict $gamma] } {
	return -code error [expectDouble $gamma]
    }
    if { $gamma <= 0.0 || $gamma >= 1.0 } {
	return -code error -errorcode [list romberg gammaTooBig] \
	    "gamma must lie in the interval (0,1)"
    }
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a >= $b } {
	return -code error "limits of integration out of order"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set onemgamma [expr { 1. - $gamma }]
    set f [list u_powerLawLower $f \
	       [expr { $gamma / $onemgamma }] \
	       [expr { 1 / $onemgamma }] \
	       $a]
	
    set limit [expr { pow( $b - $a, $onemgamma ) }]
    set result {}
    foreach v [eval [linsert $args 0 romberg $f 0 $limit]] {
	lappend result [expr { $v / $onemgamma }]
    }
    return $result

}

#----------------------------------------------------------------------
#
# u_powerLawLower --
#	Change of variable for integrating over an interval with
#	an integrable power law singularity at the upper bound.
#
# Parameters:
#	f - Function being integrated
#	gammaover1mgamma - gamma / (1 - gamma), where gamma is the power
#	oneover1mgamma - 1 / (1 - gamma), where gamma is the power
#	b - Upper limit of integration
#	u - Changed variable u == (b-x)**(1-gamma)
#
# Results:
#	Returns u**(1/1-gamma) * f(b-u**(1/1-gamma) ).
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_powerLawUpper { f gammaover1mgamma oneover1mgamma
					 b u } {
    set cmd $f
    lappend cmd [expr { $b - pow( $u, $oneover1mgamma ) }]
    set y [eval $cmd]
    return [expr { $y * pow( $u, $gammaover1mgamma ) }]
}

#----------------------------------------------------------------------
#
# math::calculus::romberg_powerLawUpper --
#	Integrate a function with an integrable power law singularity
#	at the upper bound
#
# Usage:
#	romberg_powerLawLower gamma f a b ?-option value...?
#
# Parameters:
#	gamma - Power (0<gamma<1) of the singularity
#	f - Function to integrate.  Must be a single Tcl command,
#	    to which will be appended the abscissa at which the function
#	    should be evaluated.  f is expected to have an integrable
#	    power law singularity at the upper endpoint; that is, the
#	    integrand is expected to diverge as (b-x)**gamma.
#	a - One bound of the interval
#	b - The other bound of the interval.  a and b need not be in
#	    ascending order.
#
# Options:
#	-abserror ABSERROR
#		Requests that the integration be performed to make
#		the estimated absolute error of the integral less than
#		the given value.  Default is 1.e-10.
#	-relerror RELERROR
#		Requests that the integration be performed to make
#		the estimated absolute error of the integral less than
#		the given value.  Default is 1.e-6.
#	-degree N
#		Specifies the degree of the polynomial that will be
#		used to extrapolate to a zero step size.  -degree 0
#		requests integration with the midpoint rule; -degree 1
#		is equivalent to Simpson's 3/8 rule; higher degrees
#		are difficult to describe but (within reason) give
#		faster convergence for smooth functions.  Default is
#		-degree 4.
#	-maxiter N
#		Specifies the maximum number of triplings of the
#		number of steps to take in integration.  At most
#		3**N function evaluations will be performed in
#		integrating with -maxiter N.  The integration
#		will terminate at that time, even if the result
#		satisfies neither the -relerror nor -abserror tests.
#
# Results:
#	Returns a two-element list.  The first element is the estimated
#	value of the integral; the second is the estimated absolute
#	error of the value.
#
# The 'romberg_PowerLawUpper' procedure is a wrapper for 'romberg'
# for integrating a function with an integrable power law singularity
# at the upper bound of the interval.  It works by making the change
# of variable u = (b-x)**(1-gamma).
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_powerLawUpper { gamma f a b args } {
    if { ![string is double -strict $gamma] } {
	return -code error [expectDouble $gamma]
    }
    if { $gamma <= 0.0 || $gamma >= 1.0 } {
	return -code error -errorcode [list romberg gammaTooBig] \
	    "gamma must lie in the interval (0,1)"
    }
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a >= $b } {
	return -code error "limits of integration out of order"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set onemgamma [expr { 1. - $gamma }]
    set f [list u_powerLawUpper $f \
	       [expr { $gamma / $onemgamma }] \
	       [expr { 1. / $onemgamma }] \
	       $b]
	
    set limit [expr { pow( $b - $a, $onemgamma ) }]
    set result {}
    foreach v [eval [linsert $args 0 romberg $f 0 $limit]] {
	lappend result [expr { $v / $onemgamma }]
    }
    return $result

}

#----------------------------------------------------------------------
#
# u_expUpper --
#
#	Change of variable to integrate a function that decays
#	exponentially.
#
# Parameters:
#	f - Function to integrate
#	u - Changed variable u = exp(-x)
#
# Results:
#	Returns (1/u)*f(-log(u))
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_expUpper { f u } {
    set cmd $f
    lappend cmd [expr { -log($u) }]
    set y [eval $cmd]
    return [expr { $y / $u }]
}

#----------------------------------------------------------------------
#
# romberg_expUpper --
#
#	Integrate a function that decays exponentially over a
#	half-infinite interval.
#
# Parameters:
#	Same as romberg.  The upper limit of integration, 'b',
#	is expected to be very large.
#
# Results:
#	Same as romberg.
#
# The romberg_expUpper function operates by making the change of
# variable, u = exp(-x).
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_expUpper { f a b args } {
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a >= $b } {
	return -code error "limits of integration out of order"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set f [list u_expUpper $f]
    return [eval [linsert $args 0 \
		      romberg $f [expr {exp(-$b)}] [expr {exp(-$a)}]]]
}

#----------------------------------------------------------------------
#
# u_expLower --
#
#	Change of variable to integrate a function that grows
#	exponentially.
#
# Parameters:
#	f - Function to integrate
#	u - Changed variable u = exp(x)
#
# Results:
#	Returns (1/u)*f(log(u))
#
# Side effects:
#	Whatever f does.
#
#----------------------------------------------------------------------

proc ::math::calculus::u_expLower { f u } {
    set cmd $f
    lappend cmd [expr { log($u) }]
    set y [eval $cmd]
    return [expr { $y / $u }]
}

#----------------------------------------------------------------------
#
# romberg_expLower --
#
#	Integrate a function that grows exponentially over a
#	half-infinite interval.
#
# Parameters:
#	Same as romberg.  The lower limit of integration, 'a',
#	is expected to be very large and negative.
#
# Results:
#	Same as romberg.
#
# The romberg_expUpper function operates by making the change of
# variable, u = exp(x).
#
#----------------------------------------------------------------------

proc ::math::calculus::romberg_expLower { f a b args } {
    if { ![string is double -strict $a] } {
	return -code error [expectDouble $a]
    }
    if { ![string is double -strict $b] } {
	return -code error [expectDouble $b]
    }
    if { $a >= $b } {
	return -code error "limits of integration out of order"
    }
    set f [lreplace $f 0 0 [uplevel 1 [list namespace which [lindex $f 0]]]]
    set f [list u_expLower $f]
    return [eval [linsert $args 0 \
		      romberg $f [expr {exp($a)}] [expr {exp($b)}]]]
}


# regula_falsi --
#    Compute the zero of a function via regula falsi
# Arguments:
#    f       Name of the procedure/command that evaluates the function
#    xb      Start of the interval that brackets the zero
#    xe      End of the interval that brackets the zero
#    eps     Relative error that is allowed (default: 1.0e-4)
# Result:
#    Estimate of the zero, such that the estimated (!)
#    error < eps * abs(xe-xb)
# Note:
#    f(xb)*f(xe) must be negative and eps must be positive
#
proc ::math::calculus::regula_falsi { f xb xe {eps 1.0e-4} } {
    if { $eps <= 0.0 } {
       return -code error "Relative error must be positive"
    }

    set fb [$f $xb]
    set fe [$f $xe]

    if { $fb * $fe > 0.0 } {
       return -code error "Interval must be chosen such that the \
function has a different sign at the beginning than at the end"
    }

    set max_error [expr {$eps * abs($xe-$xb)}]
    set interval  [expr {abs($xe-$xb)}]

    while { $interval > $max_error } {
       set coeff [expr {($fe-$fb)/($xe-$xb)}]
       set xi    [expr {$xb-$fb/$coeff}]
       set fi    [$f $xi]

       if { $fi == 0.0 } {
          break
       }
       set diff1 [expr {abs($xe-$xi)}]
       set diff2 [expr {abs($xb-$xi)}]
       if { $diff1 > $diff2 } {
          set interval $diff2
       } else {
          set interval $diff1
       }

       if { $fb*$fi < 0.0 } {
          set xe $xi
          set fe $fi
       } else {
          set xb $xi
          set fb $fi
       }
    }

    return $xi
}

#

# qk15_basic --
#     Apply the QK15 rule to a single interval and return all results
#
# Arguments:
#     f             Function to integrate (name of procedure)
#     xstart        Start of the interval
#     xend          End of the interval
#
# Returns:
#     List of the following:
#     result        Estimated integral (I) of function f
#     abserr        Estimate of the absolute error in "result"
#     resabs        Estimated integral of the absolute value of f
#     resasc        Estimated integral of abs(f - I/(xend-xstart))
#
# Note:
#     Translation of the 15-point Gauss-Kronrod rule (QK15) as found
#     in the SLATEC library (QUADPACK) into Tcl.
#
namespace eval ::math::calculus {
    variable qk15_xgk
    variable qk15_wgk
    variable qk15_wg

    set qk15_xgk {
           0.9914553711208126e+00    0.9491079123427585e+00
           0.8648644233597691e+00    0.7415311855993944e+00
           0.5860872354676911e+00    0.4058451513773972e+00
           0.2077849550078985e+00    0.0e+00               }
    set qk15_wgk {
           0.2293532201052922e-01    0.6309209262997855e-01
           0.1047900103222502e+00    0.1406532597155259e+00
           0.1690047266392679e+00    0.1903505780647854e+00
           0.2044329400752989e+00    0.2094821410847278e+00}
    set qk15_wg {
           0.1294849661688697e+00    0.2797053914892767e+00
           0.3818300505051189e+00    0.4179591836734694e+00}
}

if {[package vsatisfies [package present Tcl] 8.5]} {
    proc ::math::calculus::Min {a b} { expr {min ($a, $b)} }
    proc ::math::calculus::Max {a b} { expr {max ($a, $b)} }
} else {
    proc ::math::calculus::Min {a b} { if {$a < $b} { return $a } else { return $b }}
    proc ::math::calculus::Max {a b} { if {$a > $b} { return $a } else { return $b }}
}

proc ::math::calculus::qk15_basic {xstart xend func} {
    variable qk15_wg
    variable qk15_wgk
    variable qk15_xgk

    #
    # Use fixed values for epmach and uflow:
    # - epmach is the largest relative spacing.
    # - uflow is the smallest positive magnitude.

    set epmach [expr {2.3e-308}]
    set uflow  [expr {1.2e-16}]

    set centr  [expr {0.5e+00*($xstart+$xend)}]
    set hlgth  [expr {0.5e+00*($xend-$xstart)}]
    set dhlgth [expr {abs($hlgth)}]

    #
    # Compute the 15-point Kronrod approximation to
    # the integral, and estimate the absolute error.
    #
    set fc     [uplevel 2 $func $centr]
    set resg   [expr {$fc*[lindex $qk15_wg 3]}]
    set resk   [expr {$fc*[lindex $qk15_wgk 7]}]
    set resabs [expr {abs($resk)}]

    set fv1    [lrepeat 7 0.0]
    set fv2    [lrepeat 7 0.0]

    for {set j 0} {$j < 3} {incr j} {
        set jtw [expr {$j*2 +1}]
        set absc [expr {$hlgth*[lindex $qk15_xgk $jtw]}]
        set fval1 [uplevel 2 $func [expr {$centr-$absc}]]
        set fval2 [uplevel 2 $func [expr {$centr+$absc}]]
        lset fv1 $jtw $fval1
        lset fv2 $jtw $fval2
        set fsum [expr {$fval1+$fval2}]
        set resg [expr {$resg+[lindex $qk15_wg $j]*$fsum}]
        set resk [expr {$resk+[lindex $qk15_wgk $jtw]*$fsum}]
        set resabs [expr {$resabs+[lindex $qk15_wgk $jtw]*(abs($fval1)+abs($fval2))}]
    }
    for {set j 0} {$j < 4} {incr j} {
        set jtwm1 [expr {$j*2}]
        set absc [expr {$hlgth*[lindex $qk15_xgk $jtwm1]}]
        set fval1 [uplevel 2 $func [expr {$centr-$absc}]]
        set fval2 [uplevel 2 $func [expr {$centr+$absc}]]
        lset fv1 $jtwm1 $fval1
        lset fv2 $jtwm1 $fval2
        set fsum [expr {$fval1+$fval2}]
        set resk [expr {$resk+[lindex $qk15_wgk $jtwm1]*$fsum}]
        set resabs [expr {$resabs+[lindex $qk15_wgk $jtwm1]*(abs($fval1)+abs($fval2))}]
    }

    set reskh [expr {$resk*0.5e+00}]
    set resasc [expr {[lindex $qk15_wgk 7]*abs($fc-$reskh)}]

    for {set j 0} {$j < 7} {incr j} {
        set wgk    [lindex $qk15_wgk $j]
        set FV1    [lindex $fv1      $j]
        set FV2    [lindex $fv2      $j]
        set resasc [expr {$resasc+$wgk*(abs($FV1-$reskh)+abs($FV2-$reskh))}]
    }

    set result [expr {$resk*$hlgth}]
    set resabs [expr {$resabs*$dhlgth}]
    set resasc [expr {$resasc*$dhlgth}]
    set abserr [expr {abs(($resk-$resg)*$hlgth)}]
    if { $resasc != 0.0e+00 && $abserr != 0.0e+00 } {
        set abserr [expr {$resasc*[Min 0.1e+01 [expr {pow((0.2e+3*$abserr/$resasc),1.5e+00)}]]}]
    }
    if { $resabs > $uflow/(0.5e+02*$epmach) } {
        set abserr [Max [expr {($epmach*0.5e+02)*$resabs}] $abserr]
    }

    return [list $result $abserr $resabs $resasc]
}

# qk15 --
#     Apply the QK15 rule to an interval and return the estimated integral
#
# Arguments:
#     xstart        Start of the interval
#     xend          End of the interval
#     func          Function to integrate (name of procedure)
#     n             Number of subintervals (default: 1)
#
# Returns:
#     Estimated integral of function func
#
proc ::math::calculus::qk15 {xstart xend func {n 1}} {
    if { $n == 1 } {
        return [lindex [qk15_basic $xstart $xend $func] 0]
    } else {
        set dx [expr {($xend-$xstart)/double($n)}]
        set result 0.0
        for {set i 0} {$i < $n} {incr i} {
            set xb [expr {$xstart + $dx * $i}]
            set xe [expr {$xstart + $dx * ($i+1)}]

            set result [expr {$result + [lindex [qk15_basic $xb $xe $func] 0]}]
        }
    }

    return $result
}

# qk15_detailed --
#     Apply the QK15 rule to an interval and return the estimated integral
#     as well as the other values
#
# Arguments:
#     xstart        Start of the interval
#     xend          End of the interval
#     func          Function to integrate (name of procedure)
#     n             Number of subintervals (default: 1)
#
# Returns:
#     List of the following:
#     result        Estimated integral (I) of function func
#     abserr        Estimate of the absolute error in "result"
#     resabs        Estimated integral of the absolute value of f
#     resasc        Estimated integral of abs(f - I/(xend-xstart))
#
proc ::math::calculus::qk15_detailed {xstart xend func {n 1}} {
    if { $n == 1 } {
        return [qk15_basic $xstart $xend $func]
    } else {
        set dx [expr {($xend-$xstart)/double($n)}]
        set result 0.0
        set abserr 0.0
        set resabs 0.0
        set resasc 0.0
        for {set i 0} {$i < $n} {incr i} {
            set xb [expr {$xstart + $dx * $i}]
            set xe [expr {$xstart + $dx * ($i+1)}]

            foreach {dresult dabserr dresabs dresasc} [qk15_basic $xb $xe $func] break
            set result [expr {$result + $dresult}]
            set abserr [expr {$abserr + $dabserr}]
            set resabs [expr {$resabs + $dresabs}]
            set resasc [expr {$resasc + $dresasc}]
        }
    }

    return [list $result $abserr $resabs $resasc]
}