~ubuntu-branches/ubuntu/dapper/gsl-ref-html/dapper

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
1646
1647
1648
1649
1650
1651
1652
1653
1654
1655
1656
1657
1658
1659
1660
1661
1662
1663
1664
1665
1666
1667
1668
1669
1670
1671
1672
1673
1674
1675
1676
1677
1678
1679
1680
1681
1682
1683
1684
1685
1686
1687
1688
1689
1690
1691
1692
1693
1694
1695
1696
1697
1698
1699
1700
1701
1702
1703
1704
1705
1706
1707
1708
1709
1710
1711
1712
1713
1714
1715
1716
1717
1718
1719
1720
1721
1722
1723
1724
1725
1726
1727
1728
1729
1730
1731
1732
1733
1734
1735
1736
1737
1738
1739
1740
1741
1742
1743
1744
1745
1746
1747
1748
1749
1750
1751
1752
1753
1754
1755
1756
1757
1758
1759
1760
1761
1762
1763
1764
1765
1766
1767
1768
1769
1770
1771
1772
1773
1774
1775
1776
1777
1778
1779
1780
1781
1782
1783
1784
1785
1786
1787
1788
1789
1790
1791
1792
1793
1794
1795
1796
1797
1798
1799
1800
1801
1802
1803
1804
1805
1806
1807
1808
1809
1810
1811
1812
1813
1814
1815
1816
1817
1818
1819
1820
1821
1822
1823
1824
1825
1826
1827
1828
1829
1830
1831
1832
1833
1834
1835
1836
1837
1838
1839
1840
1841
1842
1843
1844
1845
1846
1847
1848
1849
1850
1851
1852
1853
1854
1855
1856
1857
1858
1859
1860
1861
1862
1863
1864
1865
1866
1867
1868
1869
1870
1871
1872
1873
1874
1875
1876
1877
1878
1879
1880
1881
1882
1883
1884
1885
1886
1887
1888
1889
1890
1891
1892
1893
1894
1895
1896
1897
1898
1899
1900
1901
1902
1903
1904
1905
1906
1907
1908
1909
1910
1911
1912
1913
1914
1915
1916
1917
1918
1919
1920
1921
1922
1923
1924
1925
1926
1927
1928
1929
1930
1931
1932
1933
1934
1935
1936
1937
1938
1939
1940
1941
1942
1943
1944
1945
1946
1947
1948
1949
1950
1951
1952
1953
1954
1955
1956
1957
1958
1959
1960
1961
1962
1963
1964
1965
1966
1967
1968
1969
1970
1971
1972
1973
1974
1975
1976
1977
1978
1979
1980
1981
1982
1983
1984
1985
1986
1987
1988
1989
1990
1991
1992
1993
1994
1995
1996
1997
1998
1999
2000
2001
2002
2003
2004
2005
2006
2007
2008
2009
2010
2011
2012
2013
2014
2015
2016
2017
2018
2019
2020
2021
2022
2023
2024
2025
2026
2027
2028
2029
2030
2031
2032
2033
2034
2035
2036
2037
2038
2039
2040
2041
2042
2043
2044
2045
2046
2047
2048
2049
2050
2051
2052
2053
2054
2055
2056
2057
2058
2059
2060
2061
2062
2063
2064
2065
2066
2067
2068
2069
2070
2071
2072
2073
2074
2075
2076
2077
2078
2079
2080
2081
2082
2083
2084
2085
2086
2087
2088
2089
2090
2091
2092
2093
2094
2095
2096
2097
2098
2099
2100
2101
2102
2103
2104
2105
2106
2107
2108
2109
2110
2111
2112
2113
2114
2115
2116
2117
2118
2119
2120
2121
2122
2123
2124
2125
2126
2127
2128
2129
2130
2131
2132
2133
2134
2135
2136
2137
2138
2139
2140
2141
2142
2143
2144
2145
2146
2147
2148
2149
2150
2151
2152
2153
2154
2155
2156
2157
2158
2159
2160
2161
2162
2163
2164
2165
2166
2167
2168
2169
2170
2171
2172
2173
2174
2175
2176
2177
2178
2179
2180
2181
2182
2183
2184
2185
2186
2187
2188
2189
2190
2191
<HTML>
<HEAD>
<!-- This HTML file has been created by texi2html 1.54+ (gsl)
     from /home/bjg/gsl.redhat/doc/gsl-ref.texi on 14 September 2005 -->

<TITLE>GNU Scientific Library -- Reference Manual - Vectors and Matrices</TITLE>
<link href="gsl-ref_9.html" rel=Next>
<link href="gsl-ref_7.html" rel=Previous>
<link href="gsl-ref_toc.html" rel=ToC>

</HEAD>
<BODY>
<p>Go to the <A HREF="gsl-ref_1.html">first</A>, <A HREF="gsl-ref_7.html">previous</A>, <A HREF="gsl-ref_9.html">next</A>, <A HREF="gsl-ref_50.html">last</A> section, <A HREF="gsl-ref_toc.html">table of contents</A>.
<P><HR><P>


<H1><A NAME="SEC155" HREF="gsl-ref_toc.html#TOC155">Vectors and Matrices</A></H1>
<P>
<A NAME="IDX809"></A>
<A NAME="IDX810"></A>
<A NAME="IDX811"></A>

</P>
<P>
The functions described in this chapter provide a simple vector and
matrix interface to ordinary C arrays. The memory management of these
arrays is implemented using a single underlying type, known as a
block. By writing your functions in terms of vectors and matrices you
can pass a single structure containing both data and dimensions as an
argument without needing additional function parameters.  The structures
are compatible with the vector and matrix formats used by BLAS
routines.

</P>



<H2><A NAME="SEC156" HREF="gsl-ref_toc.html#TOC156">Data types</A></H2>

<P>
All the functions are available for each of the standard data-types.
The versions for <CODE>double</CODE> have the prefix <CODE>gsl_block</CODE>,
<CODE>gsl_vector</CODE> and <CODE>gsl_matrix</CODE>.  Similarly the versions for
single-precision <CODE>float</CODE> arrays have the prefix
<CODE>gsl_block_float</CODE>, <CODE>gsl_vector_float</CODE> and
<CODE>gsl_matrix_float</CODE>.  The full list of available types is given
below,

</P>

<PRE>
gsl_block                       double         
gsl_block_float                 float         
gsl_block_long_double           long double   
gsl_block_int                   int           
gsl_block_uint                  unsigned int  
gsl_block_long                  long          
gsl_block_ulong                 unsigned long 
gsl_block_short                 short         
gsl_block_ushort                unsigned short
gsl_block_char                  char          
gsl_block_uchar                 unsigned char 
gsl_block_complex               complex double        
gsl_block_complex_float         complex float         
gsl_block_complex_long_double   complex long double   
</PRE>

<P>
Corresponding types exist for the <CODE>gsl_vector</CODE> and
<CODE>gsl_matrix</CODE> functions.

</P>



<H2><A NAME="SEC157" HREF="gsl-ref_toc.html#TOC157">Blocks</A></H2>

<P>
For consistency all memory is allocated through a <CODE>gsl_block</CODE>
structure.  The structure contains two components, the size of an area of
memory and a pointer to the memory.  The <CODE>gsl_block</CODE> structure looks
like this,

</P>

<PRE>
typedef struct
{
  size_t size;
  double * data;
} gsl_block;
</PRE>

<P>
Vectors and matrices are made by <I>slicing</I> an underlying block. A
slice is a set of elements formed from an initial offset and a
combination of indices and step-sizes. In the case of a matrix the
step-size for the column index represents the row-length.  The step-size
for a vector is known as the <I>stride</I>.

</P>
<P>
The functions for allocating and deallocating blocks are defined in
<TT>`gsl_block.h'</TT>

</P>



<H3><A NAME="SEC158" HREF="gsl-ref_toc.html#TOC158">Block allocation</A></H3>

<P>
The functions for allocating memory to a block follow the style of
<CODE>malloc</CODE> and <CODE>free</CODE>.  In addition they also perform their own
error checking.  If there is insufficient memory available to allocate a
block then the functions call the GSL error handler (with an error
number of <CODE>GSL_ENOMEM</CODE>) in addition to returning a null
pointer.  Thus if you use the library error handler to abort your program
then it isn't necessary to check every <CODE>alloc</CODE>.

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_block * <B>gsl_block_alloc</B> <I>(size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX812"></A>
This function allocates memory for a block of <VAR>n</VAR> double-precision
elements, returning a pointer to the block struct.  The block is not
initialized and so the values of its elements are undefined.  Use the
function <CODE>gsl_block_calloc</CODE> if you want to ensure that all the
elements are initialized to zero.

</P>
<P>
A null pointer is returned if insufficient memory is available to create
the block.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_block * <B>gsl_block_calloc</B> <I>(size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX813"></A>
This function allocates memory for a block and initializes all the
elements of the block to zero.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_block_free</B> <I>(gsl_block * <VAR>b</VAR>)</I>
<DD><A NAME="IDX814"></A>
This function frees the memory used by a block <VAR>b</VAR> previously
allocated with <CODE>gsl_block_alloc</CODE> or <CODE>gsl_block_calloc</CODE>.
</DL>

</P>


<H3><A NAME="SEC159" HREF="gsl-ref_toc.html#TOC159">Reading and writing blocks</A></H3>

<P>
The library provides functions for reading and writing blocks to a file
as binary data or formatted text.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_block_fwrite</B> <I>(FILE * <VAR>stream</VAR>, const gsl_block * <VAR>b</VAR>)</I>
<DD><A NAME="IDX815"></A>
This function writes the elements of the block <VAR>b</VAR> to the stream
<VAR>stream</VAR> in binary format.  The return value is 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem writing to the file.  Since the
data is written in the native binary format it may not be portable
between different architectures.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_block_fread</B> <I>(FILE * <VAR>stream</VAR>, gsl_block * <VAR>b</VAR>)</I>
<DD><A NAME="IDX816"></A>
This function reads into the block <VAR>b</VAR> from the open stream
<VAR>stream</VAR> in binary format.  The block <VAR>b</VAR> must be preallocated
with the correct length since the function uses the size of <VAR>b</VAR> to
determine how many bytes to read.  The return value is 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem reading from the file.  The
data is assumed to have been written in the native binary format on the
same architecture.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_block_fprintf</B> <I>(FILE * <VAR>stream</VAR>, const gsl_block * <VAR>b</VAR>, const char * <VAR>format</VAR>)</I>
<DD><A NAME="IDX817"></A>
This function writes the elements of the block <VAR>b</VAR> line-by-line to
the stream <VAR>stream</VAR> using the format specifier <VAR>format</VAR>, which
should be one of the <CODE>%g</CODE>, <CODE>%e</CODE> or <CODE>%f</CODE> formats for
floating point numbers and <CODE>%d</CODE> for integers.  The function returns
0 for success and <CODE>GSL_EFAILED</CODE> if there was a problem writing to
the file.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_block_fscanf</B> <I>(FILE * <VAR>stream</VAR>, gsl_block * <VAR>b</VAR>)</I>
<DD><A NAME="IDX818"></A>
This function reads formatted data from the stream <VAR>stream</VAR> into the
block <VAR>b</VAR>.  The block <VAR>b</VAR> must be preallocated with the correct
length since the function uses the size of <VAR>b</VAR> to determine how many
numbers to read.  The function returns 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem reading from the file.
</DL>

</P>


<H3><A NAME="SEC160" HREF="gsl-ref_toc.html#TOC160">Example programs for blocks</A></H3>

<P>
The following program shows how to allocate a block,

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_block.h&#62;

int
main (void)
{
  gsl_block * b = gsl_block_alloc (100);
  
  printf ("length of block = %u\n", b-&#62;size);
  printf ("block data address = %#x\n", b-&#62;data);

  gsl_block_free (b);
  return 0;
}
</PRE>

<P>
Here is the output from the program,

</P>

<PRE>
length of block = 100
block data address = 0x804b0d8
</PRE>



<H2><A NAME="SEC161" HREF="gsl-ref_toc.html#TOC161">Vectors</A></H2>
<P>
<A NAME="IDX819"></A>
<A NAME="IDX820"></A>

</P>
<P>
Vectors are defined by a <CODE>gsl_vector</CODE> structure which describes a
slice of a block.  Different vectors can be created which point to the
same block.  A vector slice is a set of equally-spaced elements of an
area of memory.

</P>
<P>
The <CODE>gsl_vector</CODE> structure contains five components, the
<I>size</I>, the <I>stride</I>, a pointer to the memory where the elements
are stored, <VAR>data</VAR>, a pointer to the block owned by the vector,
<VAR>block</VAR>, if any, and an ownership flag, <VAR>owner</VAR>.  The structure
is very simple and looks like this,

</P>

<PRE>
typedef struct
{
  size_t size;
  size_t stride;
  double * data;
  gsl_block * block;
  int owner;
} gsl_vector;
</PRE>

<P>
The <VAR>size</VAR> is simply the number of vector elements.  The range of
valid indices runs from 0 to <CODE>size-1</CODE>.  The <VAR>stride</VAR> is the
step-size from one element to the next in physical memory, measured in
units of the appropriate datatype.  The pointer <VAR>data</VAR> gives the
location of the first element of the vector in memory.  The pointer
<VAR>block</VAR> stores the location of the memory block in which the vector
elements are located (if any).  If the vector owns this block then the
<VAR>owner</VAR> field is set to one and the block will be deallocated when the
vector is freed.  If the vector points to a block owned by another
object then the <VAR>owner</VAR> field is zero and any underlying block will not be
deallocated with the vector.

</P>
<P>
The functions for allocating and accessing vectors are defined in
<TT>`gsl_vector.h'</TT>

</P>



<H3><A NAME="SEC162" HREF="gsl-ref_toc.html#TOC162">Vector allocation</A></H3>

<P>
The functions for allocating memory to a vector follow the style of
<CODE>malloc</CODE> and <CODE>free</CODE>.  In addition they also perform their own
error checking.  If there is insufficient memory available to allocate a
vector then the functions call the GSL error handler (with an error
number of <CODE>GSL_ENOMEM</CODE>) in addition to returning a null
pointer.  Thus if you use the library error handler to abort your program
then it isn't necessary to check every <CODE>alloc</CODE>.

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector * <B>gsl_vector_alloc</B> <I>(size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX821"></A>
This function creates a vector of length <VAR>n</VAR>, returning a pointer to
a newly initialized vector struct. A new block is allocated for the
elements of the vector, and stored in the <VAR>block</VAR> component of the
vector struct.  The block is "owned" by the vector, and will be
deallocated when the vector is deallocated.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector * <B>gsl_vector_calloc</B> <I>(size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX822"></A>
This function allocates memory for a vector of length <VAR>n</VAR> and
initializes all the elements of the vector to zero.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_vector_free</B> <I>(gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX823"></A>
This function frees a previously allocated vector <VAR>v</VAR>.  If the
vector was created using <CODE>gsl_vector_alloc</CODE> then the block
underlying the vector will also be deallocated.  If the vector has been
created from another object then the memory is still owned by that
object and will not be deallocated.
</DL>

</P>


<H3><A NAME="SEC163" HREF="gsl-ref_toc.html#TOC163">Accessing vector elements</A></H3>
<P>
<A NAME="IDX824"></A>
<A NAME="IDX825"></A>
<A NAME="IDX826"></A>
<A NAME="IDX827"></A>
<A NAME="IDX828"></A>

</P>
<P>
Unlike FORTRAN compilers, C compilers do not usually provide
support for range checking of vectors and matrices.  Range checking is
available in the GNU C Compiler bounds-checking extension, but it is not
part of the default installation of GCC.  The functions
<CODE>gsl_vector_get</CODE> and <CODE>gsl_vector_set</CODE> can perform portable
range checking for you and report an error if you attempt to access
elements outside the allowed range.

</P>
<P>
The functions for accessing the elements of a vector or matrix are
defined in <TT>`gsl_vector.h'</TT> and declared <CODE>extern inline</CODE> to
eliminate function-call overhead.  You must compile your program with
the macro <CODE>HAVE_INLINE</CODE> defined to use these functions.  

</P>
<P>
If necessary you can turn off range checking completely without
modifying any source files by recompiling your program with the
preprocessor definition <CODE>GSL_RANGE_CHECK_OFF</CODE>.  Provided your
compiler supports inline functions the effect of turning off range
checking is to replace calls to <CODE>gsl_vector_get(v,i)</CODE> by
<CODE>v-&#62;data[i*v-&#62;stride]</CODE> and calls to <CODE>gsl_vector_set(v,i,x)</CODE> by
<CODE>v-&#62;data[i*v-&#62;stride]=x</CODE>.  Thus there should be no performance
penalty for using the range checking functions when range checking is
turned off.

</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_vector_get</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX829"></A>
This function returns the <VAR>i</VAR>-th element of a vector <VAR>v</VAR>.  If
<VAR>i</VAR> lies outside the allowed range of 0 to <VAR>n</VAR>-1 then the error
handler is invoked and 0 is returned.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_vector_set</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>i</VAR>, double <VAR>x</VAR>)</I>
<DD><A NAME="IDX830"></A>
This function sets the value of the <VAR>i</VAR>-th element of a vector
<VAR>v</VAR> to <VAR>x</VAR>.  If <VAR>i</VAR> lies outside the allowed range of 0 to
<VAR>n</VAR>-1 then the error handler is invoked.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> double * <B>gsl_vector_ptr</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX831"></A>
<DT><U>Function:</U> const double * <B>gsl_vector_const_ptr</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX832"></A>
These functions return a pointer to the <VAR>i</VAR>-th element of a vector
<VAR>v</VAR>.  If <VAR>i</VAR> lies outside the allowed range of 0 to <VAR>n</VAR>-1
then the error handler is invoked and a null pointer is returned.
</DL>

</P>


<H3><A NAME="SEC164" HREF="gsl-ref_toc.html#TOC164">Initializing vector elements</A></H3>
<P>
<A NAME="IDX833"></A>
<A NAME="IDX834"></A>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_vector_set_all</B> <I>(gsl_vector * <VAR>v</VAR>, double <VAR>x</VAR>)</I>
<DD><A NAME="IDX835"></A>
This function sets all the elements of the vector <VAR>v</VAR> to the value
<VAR>x</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_vector_set_zero</B> <I>(gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX836"></A>
This function sets all the elements of the vector <VAR>v</VAR> to zero.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_set_basis</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX837"></A>
This function makes a basis vector by setting all the elements of the
vector <VAR>v</VAR> to zero except for the <VAR>i</VAR>-th element which is set to
one.
</DL>

</P>


<H3><A NAME="SEC165" HREF="gsl-ref_toc.html#TOC165">Reading and writing vectors</A></H3>

<P>
The library provides functions for reading and writing vectors to a file
as binary data or formatted text.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_fwrite</B> <I>(FILE * <VAR>stream</VAR>, const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX838"></A>
This function writes the elements of the vector <VAR>v</VAR> to the stream
<VAR>stream</VAR> in binary format.  The return value is 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem writing to the file.  Since the
data is written in the native binary format it may not be portable
between different architectures.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_fread</B> <I>(FILE * <VAR>stream</VAR>, gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX839"></A>
This function reads into the vector <VAR>v</VAR> from the open stream
<VAR>stream</VAR> in binary format.  The vector <VAR>v</VAR> must be preallocated
with the correct length since the function uses the size of <VAR>v</VAR> to
determine how many bytes to read.  The return value is 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem reading from the file.  The
data is assumed to have been written in the native binary format on the
same architecture.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_fprintf</B> <I>(FILE * <VAR>stream</VAR>, const gsl_vector * <VAR>v</VAR>, const char * <VAR>format</VAR>)</I>
<DD><A NAME="IDX840"></A>
This function writes the elements of the vector <VAR>v</VAR> line-by-line to
the stream <VAR>stream</VAR> using the format specifier <VAR>format</VAR>, which
should be one of the <CODE>%g</CODE>, <CODE>%e</CODE> or <CODE>%f</CODE> formats for
floating point numbers and <CODE>%d</CODE> for integers.  The function returns
0 for success and <CODE>GSL_EFAILED</CODE> if there was a problem writing to
the file.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_fscanf</B> <I>(FILE * <VAR>stream</VAR>, gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX841"></A>
This function reads formatted data from the stream <VAR>stream</VAR> into the
vector <VAR>v</VAR>.  The vector <VAR>v</VAR> must be preallocated with the correct
length since the function uses the size of <VAR>v</VAR> to determine how many
numbers to read.  The function returns 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem reading from the file.
</DL>

</P>


<H3><A NAME="SEC166" HREF="gsl-ref_toc.html#TOC166">Vector views</A></H3>

<P>
In addition to creating vectors from slices of blocks it is also
possible to slice vectors and create vector views.  For example, a
subvector of another vector can be described with a view, or two views
can be made which provide access to the even and odd elements of a
vector.

</P>
<P>
A vector view is a temporary object, stored on the stack, which can be
used to operate on a subset of vector elements.  Vector views can be
defined for both constant and non-constant vectors, using separate types
that preserve constness.  A vector view has the type
<CODE>gsl_vector_view</CODE> and a constant vector view has the type
<CODE>gsl_vector_const_view</CODE>.  In both cases the elements of the view
can be accessed as a <CODE>gsl_vector</CODE> using the <CODE>vector</CODE> component
of the view object.  A pointer to a vector of type <CODE>gsl_vector *</CODE>
or <CODE>const gsl_vector *</CODE> can be obtained by taking the address of
this component with the <CODE>&#38;</CODE> operator.  

</P>
<P>
When using this pointer it is important to ensure that the view itself
remains in scope--the simplest way to do so is by always writing the
pointer as <CODE>&#38;</CODE><VAR>view</VAR><CODE>.vector</CODE>, and never storing this value
in another variable.  

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_vector_subvector</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>offset</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX842"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_vector_const_subvector</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t <VAR>offset</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX843"></A>
These functions return a vector view of a subvector of another vector
<VAR>v</VAR>.  The start of the new vector is offset by <VAR>offset</VAR> elements
from the start of the original vector.  The new vector has <VAR>n</VAR>
elements.  Mathematically, the <VAR>i</VAR>-th element of the new vector
<VAR>v'</VAR> is given by,

</P>

<PRE>
v'(i) = v-&#62;data[(offset + i)*v-&#62;stride]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n-1</CODE>.

</P>
<P>
The <CODE>data</CODE> pointer of the returned vector struct is set to null if
the combined parameters (<VAR>offset</VAR>,<VAR>n</VAR>) overrun the end of the
original vector.

</P>
<P>
The new vector is only a view of the block underlying the original
vector, <VAR>v</VAR>.  The block containing the elements of <VAR>v</VAR> is not
owned by the new vector.  When the view goes out of scope the original
vector <VAR>v</VAR> and its block will continue to exist.  The original
memory can only be deallocated by freeing the original vector.  Of
course, the original vector should not be deallocated while the view is
still in use.

</P>
<P>
The function <CODE>gsl_vector_const_subvector</CODE> is equivalent to
<CODE>gsl_vector_subvector</CODE> but can be used for vectors which are
declared <CODE>const</CODE>.
</DL>

</P>

<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_vector_subvector_with_stride</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>offset</VAR>, size_t <VAR>stride</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX844"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_vector_const_subvector_with_stride</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t <VAR>offset</VAR>, size_t <VAR>stride</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX845"></A>
These functions return a vector view of a subvector of another vector
<VAR>v</VAR> with an additional stride argument. The subvector is formed in
the same way as for <CODE>gsl_vector_subvector</CODE> but the new vector has
<VAR>n</VAR> elements with a step-size of <VAR>stride</VAR> from one element to
the next in the original vector.  Mathematically, the <VAR>i</VAR>-th element
of the new vector <VAR>v'</VAR> is given by,

</P>

<PRE>
v'(i) = v-&#62;data[(offset + i*stride)*v-&#62;stride]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n-1</CODE>.

</P>
<P>
Note that subvector views give direct access to the underlying elements
of the original vector. For example, the following code will zero the
even elements of the vector <CODE>v</CODE> of length <CODE>n</CODE>, while leaving the
odd elements untouched,

</P>

<PRE>
gsl_vector_view v_even 
  = gsl_vector_subvector_with_stride (v, 0, 2, n/2);
gsl_vector_set_zero (&#38;v_even.vector);
</PRE>

<P>
A vector view can be passed to any subroutine which takes a vector
argument just as a directly allocated vector would be, using
<CODE>&#38;</CODE><VAR>view</VAR><CODE>.vector</CODE>.  For example, the following code
computes the norm of the odd elements of <CODE>v</CODE> using the BLAS
routine DNRM2,

</P>

<PRE>
gsl_vector_view v_odd 
  = gsl_vector_subvector_with_stride (v, 1, 2, n/2);
double r = gsl_blas_dnrm2 (&#38;v_odd.vector);
</PRE>

<P>
The function <CODE>gsl_vector_const_subvector_with_stride</CODE> is equivalent
to <CODE>gsl_vector_subvector_with_stride</CODE> but can be used for vectors
which are declared <CODE>const</CODE>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_vector_complex_real</B> <I>(gsl_vector_complex * <VAR>v</VAR>)</I>
<DD><A NAME="IDX846"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_vector_complex_const_real</B> <I>(const gsl_vector_complex * <VAR>v</VAR>)</I>
<DD><A NAME="IDX847"></A>
These functions return a vector view of the real parts of the complex
vector <VAR>v</VAR>.

</P>
<P>
The function <CODE>gsl_vector_complex_const_real</CODE> is equivalent to
<CODE>gsl_vector_complex_real</CODE> but can be used for vectors which are
declared <CODE>const</CODE>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_vector_complex_imag</B> <I>(gsl_vector_complex * <VAR>v</VAR>)</I>
<DD><A NAME="IDX848"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_vector_complex_const_imag</B> <I>(const gsl_vector_complex * <VAR>v</VAR>)</I>
<DD><A NAME="IDX849"></A>
These functions return a vector view of the imaginary parts of the
complex vector <VAR>v</VAR>.

</P>
<P>
The function <CODE>gsl_vector_complex_const_imag</CODE> is equivalent to
<CODE>gsl_vector_complex_imag</CODE> but can be used for vectors which are
declared <CODE>const</CODE>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_vector_view_array</B> <I>(double * <VAR>base</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX850"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_vector_const_view_array</B> <I>(const double * <VAR>base</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX851"></A>
These functions return a vector view of an array.  The start of the new
vector is given by <VAR>base</VAR> and has <VAR>n</VAR> elements.  Mathematically,
the <VAR>i</VAR>-th element of the new vector <VAR>v'</VAR> is given by,

</P>

<PRE>
v'(i) = base[i]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n-1</CODE>.

</P>
<P>
The array containing the elements of <VAR>v</VAR> is not owned by the new
vector view.  When the view goes out of scope the original array will
continue to exist.  The original memory can only be deallocated by
freeing the original pointer <VAR>base</VAR>.  Of course, the original array
should not be deallocated while the view is still in use.

</P>
<P>
The function <CODE>gsl_vector_const_view_array</CODE> is equivalent to
<CODE>gsl_vector_view_array</CODE> but can be used for arrays which are
declared <CODE>const</CODE>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_vector_view_array_with_stride</B> <I>(double * <VAR>base</VAR>, size_t <VAR>stride</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX852"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_vector_const_view_array_with_stride</B> <I>(const double * <VAR>base</VAR>, size_t <VAR>stride</VAR>, size_t <VAR>n</VAR>)</I>
<DD><A NAME="IDX853"></A>
These functions return a vector view of an array <VAR>base</VAR> with an
additional stride argument. The subvector is formed in the same way as
for <CODE>gsl_vector_view_array</CODE> but the new vector has <VAR>n</VAR> elements
with a step-size of <VAR>stride</VAR> from one element to the next in the
original array.  Mathematically, the <VAR>i</VAR>-th element of the new
vector <VAR>v'</VAR> is given by,

</P>

<PRE>
v'(i) = base[i*stride]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n-1</CODE>.

</P>
<P>
Note that the view gives direct access to the underlying elements of the
original array.  A vector view can be passed to any subroutine which
takes a vector argument just as a directly allocated vector would be,
using <CODE>&#38;</CODE><VAR>view</VAR><CODE>.vector</CODE>.

</P>
<P>
The function <CODE>gsl_vector_const_view_array_with_stride</CODE> is
equivalent to <CODE>gsl_vector_view_array_with_stride</CODE> but can be used
for arrays which are declared <CODE>const</CODE>.
</DL>

</P>



<H3><A NAME="SEC167" HREF="gsl-ref_toc.html#TOC167">Copying vectors</A></H3>

<P>
Common operations on vectors such as addition and multiplication are
available in the BLAS part of the library (see section <A HREF="gsl-ref_12.html#SEC215">BLAS Support</A>).  However, it is useful to have a small number of utility
functions which do not require the full BLAS code.  The following
functions fall into this category.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_memcpy</B> <I>(gsl_vector * <VAR>dest</VAR>, const gsl_vector * <VAR>src</VAR>)</I>
<DD><A NAME="IDX854"></A>
This function copies the elements of the vector <VAR>src</VAR> into the
vector <VAR>dest</VAR>.  The two vectors must have the same length.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_swap</B> <I>(gsl_vector * <VAR>v</VAR>, gsl_vector * <VAR>w</VAR>)</I>
<DD><A NAME="IDX855"></A>
This function exchanges the elements of the vectors <VAR>v</VAR> and <VAR>w</VAR>
by copying.  The two vectors must have the same length.
</DL>

</P>



<H3><A NAME="SEC168" HREF="gsl-ref_toc.html#TOC168">Exchanging elements</A></H3>

<P>
The following function can be used to exchange, or permute, the elements
of a vector.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_swap_elements</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX856"></A>
This function exchanges the <VAR>i</VAR>-th and <VAR>j</VAR>-th elements of the
vector <VAR>v</VAR> in-place.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_reverse</B> <I>(gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX857"></A>
This function reverses the order of the elements of the vector <VAR>v</VAR>.
</DL>

</P>



<H3><A NAME="SEC169" HREF="gsl-ref_toc.html#TOC169">Vector operations</A></H3>

<P>
The following operations are only defined for real vectors.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_add</B> <I>(gsl_vector * <VAR>a</VAR>, const gsl_vector * <VAR>b</VAR>)</I>
<DD><A NAME="IDX858"></A>
This function adds the elements of vector <VAR>b</VAR> to the elements of
vector <VAR>a</VAR>, a'_i = a_i + b_i. The two vectors must have the
same length.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_sub</B> <I>(gsl_vector * <VAR>a</VAR>, const gsl_vector * <VAR>b</VAR>)</I>
<DD><A NAME="IDX859"></A>
This function subtracts the elements of vector <VAR>b</VAR> from the elements of
vector <VAR>a</VAR>, a'_i = a_i - b_i. The two vectors must have the
same length.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_mul</B> <I>(gsl_vector * <VAR>a</VAR>, const gsl_vector * <VAR>b</VAR>)</I>
<DD><A NAME="IDX860"></A>
This function multiplies the elements of vector <VAR>a</VAR> by the elements of
vector <VAR>b</VAR>, a'_i = a_i * b_i. The two vectors must have the
same length.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_div</B> <I>(gsl_vector * <VAR>a</VAR>, const gsl_vector * <VAR>b</VAR>)</I>
<DD><A NAME="IDX861"></A>
This function divides the elements of vector <VAR>a</VAR> by the elements of
vector <VAR>b</VAR>, a'_i = a_i / b_i. The two vectors must have the
same length.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_scale</B> <I>(gsl_vector * <VAR>a</VAR>, const double <VAR>x</VAR>)</I>
<DD><A NAME="IDX862"></A>
This function multiplies the elements of vector <VAR>a</VAR> by the constant
factor <VAR>x</VAR>, a'_i = x a_i.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_add_constant</B> <I>(gsl_vector * <VAR>a</VAR>, const double <VAR>x</VAR>)</I>
<DD><A NAME="IDX863"></A>
This function adds the constant value <VAR>x</VAR> to the elements of the
vector <VAR>a</VAR>, a'_i = a_i + x.
</DL>

</P>


<H3><A NAME="SEC170" HREF="gsl-ref_toc.html#TOC170">Finding maximum and minimum elements of vectors</A></H3>

<P>
<DL>
<DT><U>Function:</U> double <B>gsl_vector_max</B> <I>(const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX864"></A>
This function returns the maximum value in the vector <VAR>v</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_vector_min</B> <I>(const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX865"></A>
This function returns the minimum value in the vector <VAR>v</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_vector_minmax</B> <I>(const gsl_vector * <VAR>v</VAR>, double * <VAR>min_out</VAR>, double * <VAR>max_out</VAR>)</I>
<DD><A NAME="IDX866"></A>
This function returns the minimum and maximum values in the vector
<VAR>v</VAR>, storing them in <VAR>min_out</VAR> and <VAR>max_out</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> size_t <B>gsl_vector_max_index</B> <I>(const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX867"></A>
This function returns the index of the maximum value in the vector <VAR>v</VAR>.
When there are several equal maximum elements then the lowest index is
returned.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> size_t <B>gsl_vector_min_index</B> <I>(const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX868"></A>
This function returns the index of the minimum value in the vector <VAR>v</VAR>.
When there are several equal minimum elements then the lowest index is
returned.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_vector_minmax_index</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t * <VAR>imin</VAR>, size_t * <VAR>imax</VAR>)</I>
<DD><A NAME="IDX869"></A>
This function returns the indices of the minimum and maximum values in
the vector <VAR>v</VAR>, storing them in <VAR>imin</VAR> and <VAR>imax</VAR>. When
there are several equal minimum or maximum elements then the lowest
indices are returned.
</DL>

</P>


<H3><A NAME="SEC171" HREF="gsl-ref_toc.html#TOC171">Vector properties</A></H3>

<P>
<DL>
<DT><U>Function:</U> int <B>gsl_vector_isnull</B> <I>(const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX870"></A>
This function returns 1 if all the elements of the vector <VAR>v</VAR> are
zero, and 0 otherwise.
</DL>

</P>


<H3><A NAME="SEC172" HREF="gsl-ref_toc.html#TOC172">Example programs for vectors</A></H3>

<P>
This program shows how to allocate, initialize and read from a vector
using the functions <CODE>gsl_vector_alloc</CODE>, <CODE>gsl_vector_set</CODE> and
<CODE>gsl_vector_get</CODE>.

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_vector.h&#62;

int
main (void)
{
  int i;
  gsl_vector * v = gsl_vector_alloc (3);
  
  for (i = 0; i &#60; 3; i++)
    {
      gsl_vector_set (v, i, 1.23 + i);
    }
  
  for (i = 0; i &#60; 100; i++)
    {
      printf ("v_%d = %g\n", i, gsl_vector_get (v, i));
    }

  return 0;
}
</PRE>

<P>
Here is the output from the program.  The final loop attempts to read
outside the range of the vector <CODE>v</CODE>, and the error is trapped by
the range-checking code in <CODE>gsl_vector_get</CODE>.

</P>

<PRE>
$ ./a.out
v_0 = 1.23
v_1 = 2.23
v_2 = 3.23
gsl: vector_source.c:12: ERROR: index out of range
Default GSL error handler invoked.
Aborted (core dumped)
</PRE>

<P>
The next program shows how to write a vector to a file.

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_vector.h&#62;

int
main (void)
{
  int i; 
  gsl_vector * v = gsl_vector_alloc (100);
  
  for (i = 0; i &#60; 100; i++)
    {
      gsl_vector_set (v, i, 1.23 + i);
    }

  {  
     FILE * f = fopen ("test.dat", "w");
     gsl_vector_fprintf (f, v, "%.5g");
     fclose (f);
  }
  return 0;
}
</PRE>

<P>
After running this program the file <TT>`test.dat'</TT> should contain the
elements of <CODE>v</CODE>, written using the format specifier
<CODE>%.5g</CODE>.  The vector could then be read back in using the function
<CODE>gsl_vector_fscanf (f, v)</CODE> as follows:

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_vector.h&#62;

int
main (void)
{
  int i; 
  gsl_vector * v = gsl_vector_alloc (10);

  {  
     FILE * f = fopen ("test.dat", "r");
     gsl_vector_fscanf (f, v);
     fclose (f);
  }

  for (i = 0; i &#60; 10; i++)
    {
      printf ("%g\n", gsl_vector_get(v, i));
    }

  return 0;
}
</PRE>



<H2><A NAME="SEC173" HREF="gsl-ref_toc.html#TOC173">Matrices</A></H2>
<P>
<A NAME="IDX871"></A>
<A NAME="IDX872"></A>
<A NAME="IDX873"></A>
<A NAME="IDX874"></A>

</P>
<P>
Matrices are defined by a <CODE>gsl_matrix</CODE> structure which describes a
generalized slice of a block.  Like a vector it represents a set of
elements in an area of memory, but uses two indices instead of one.

</P>
<P>
The <CODE>gsl_matrix</CODE> structure contains six components, the two
dimensions of the matrix, a physical dimension, a pointer to the memory
where the elements of the matrix are stored, <VAR>data</VAR>, a pointer to
the block owned by the matrix <VAR>block</VAR>, if any, and an ownership
flag, <VAR>owner</VAR>.  The physical dimension determines the memory layout
and can differ from the matrix dimension to allow the use of
submatrices.  The <CODE>gsl_matrix</CODE> structure is very simple and looks
like this,

</P>

<PRE>
typedef struct
{
  size_t size1;
  size_t size2;
  size_t tda;
  double * data;
  gsl_block * block;
  int owner;
} gsl_matrix;
</PRE>

<P>
Matrices are stored in row-major order, meaning that each row of
elements forms a contiguous block in memory.  This is the standard
"C-language ordering" of two-dimensional arrays. Note that FORTRAN
stores arrays in column-major order. The number of rows is <VAR>size1</VAR>.
The range of valid row indices runs from 0 to <CODE>size1-1</CODE>.  Similarly
<VAR>size2</VAR> is the number of columns.  The range of valid column indices
runs from 0 to <CODE>size2-1</CODE>.  The physical row dimension <VAR>tda</VAR>, or
<I>trailing dimension</I>, specifies the size of a row of the matrix as
laid out in memory.

</P>
<P>
For example, in the following matrix <VAR>size1</VAR> is 3, <VAR>size2</VAR> is 4,
and <VAR>tda</VAR> is 8.  The physical memory layout of the matrix begins in
the top left hand-corner and proceeds from left to right along each row
in turn.

</P>

<PRE>
00 01 02 03 XX XX XX XX
10 11 12 13 XX XX XX XX
20 21 22 23 XX XX XX XX
</PRE>

<P>
Each unused memory location is represented by "<CODE>XX</CODE>".  The
pointer <VAR>data</VAR> gives the location of the first element of the matrix
in memory.  The pointer <VAR>block</VAR> stores the location of the memory
block in which the elements of the matrix are located (if any).  If the
matrix owns this block then the <VAR>owner</VAR> field is set to one and the
block will be deallocated when the matrix is freed.  If the matrix is
only a slice of a block owned by another object then the <VAR>owner</VAR> field is
zero and any underlying block will not be freed.

</P>
<P>
The functions for allocating and accessing matrices are defined in
<TT>`gsl_matrix.h'</TT>

</P>



<H3><A NAME="SEC174" HREF="gsl-ref_toc.html#TOC174">Matrix allocation</A></H3>

<P>
The functions for allocating memory to a matrix follow the style of
<CODE>malloc</CODE> and <CODE>free</CODE>.  They also perform their own error
checking.  If there is insufficient memory available to allocate a vector
then the functions call the GSL error handler (with an error number of
<CODE>GSL_ENOMEM</CODE>) in addition to returning a null pointer.  Thus if you
use the library error handler to abort your program then it isn't
necessary to check every <CODE>alloc</CODE>.

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_matrix * <B>gsl_matrix_alloc</B> <I>(size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX875"></A>
This function creates a matrix of size <VAR>n1</VAR> rows by <VAR>n2</VAR>
columns, returning a pointer to a newly initialized matrix struct. A new
block is allocated for the elements of the matrix, and stored in the
<VAR>block</VAR> component of the matrix struct.  The block is "owned" by the
matrix, and will be deallocated when the matrix is deallocated.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_matrix * <B>gsl_matrix_calloc</B> <I>(size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX876"></A>
This function allocates memory for a matrix of size <VAR>n1</VAR> rows by
<VAR>n2</VAR> columns and initializes all the elements of the matrix to zero.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_free</B> <I>(gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX877"></A>
This function frees a previously allocated matrix <VAR>m</VAR>.  If the
matrix was created using <CODE>gsl_matrix_alloc</CODE> then the block
underlying the matrix will also be deallocated.  If the matrix has been
created from another object then the memory is still owned by that
object and will not be deallocated.
</DL>

</P>


<H3><A NAME="SEC175" HREF="gsl-ref_toc.html#TOC175">Accessing matrix elements</A></H3>
<P>
<A NAME="IDX878"></A>
<A NAME="IDX879"></A>

</P>
<P>
The functions for accessing the elements of a matrix use the same range
checking system as vectors.  You can turn off range checking by recompiling
your program with the preprocessor definition
<CODE>GSL_RANGE_CHECK_OFF</CODE>.

</P>
<P>
The elements of the matrix are stored in "C-order", where the second
index moves continuously through memory.  More precisely, the element
accessed by the function <CODE>gsl_matrix_get(m,i,j)</CODE> and
<CODE>gsl_matrix_set(m,i,j,x)</CODE> is 

</P>

<PRE>
m-&#62;data[i * m-&#62;tda + j]
</PRE>

<P>
where <VAR>tda</VAR> is the physical row-length of the matrix.

</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_matrix_get</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX880"></A>
This function returns the (i,j)-th element of a matrix
<VAR>m</VAR>.  If <VAR>i</VAR> or <VAR>j</VAR> lie outside the allowed range of 0 to
<VAR>n1</VAR>-1 and 0 to <VAR>n2</VAR>-1 then the error handler is invoked and 0
is returned.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_set</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>, double <VAR>x</VAR>)</I>
<DD><A NAME="IDX881"></A>
This function sets the value of the (i,j)-th element of a
matrix <VAR>m</VAR> to <VAR>x</VAR>.  If <VAR>i</VAR> or <VAR>j</VAR> lies outside the
allowed range of 0 to <VAR>n1</VAR>-1 and 0 to <VAR>n2</VAR>-1 then the error
handler is invoked.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> double * <B>gsl_matrix_ptr</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX882"></A>
<DT><U>Function:</U> const double * <B>gsl_matrix_const_ptr</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX883"></A>
These functions return a pointer to the (i,j)-th element of a
matrix <VAR>m</VAR>.  If <VAR>i</VAR> or <VAR>j</VAR> lie outside the allowed range of
0 to <VAR>n1</VAR>-1 and 0 to <VAR>n2</VAR>-1 then the error handler is invoked
and a null pointer is returned.
</DL>

</P>


<H3><A NAME="SEC176" HREF="gsl-ref_toc.html#TOC176">Initializing matrix elements</A></H3>
<P>
<A NAME="IDX884"></A>
<A NAME="IDX885"></A>
<A NAME="IDX886"></A>
<A NAME="IDX887"></A>
<A NAME="IDX888"></A>
<A NAME="IDX889"></A>
<A NAME="IDX890"></A>
<A NAME="IDX891"></A>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_set_all</B> <I>(gsl_matrix * <VAR>m</VAR>, double <VAR>x</VAR>)</I>
<DD><A NAME="IDX892"></A>
This function sets all the elements of the matrix <VAR>m</VAR> to the value
<VAR>x</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_set_zero</B> <I>(gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX893"></A>
This function sets all the elements of the matrix <VAR>m</VAR> to zero.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_set_identity</B> <I>(gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX894"></A>
This function sets the elements of the matrix <VAR>m</VAR> to the
corresponding elements of the identity matrix, m(i,j) =
\delta(i,j), i.e. a unit diagonal with all off-diagonal elements zero.
This applies to both square and rectangular matrices.
</DL>

</P>


<H3><A NAME="SEC177" HREF="gsl-ref_toc.html#TOC177">Reading and writing matrices</A></H3>

<P>
The library provides functions for reading and writing matrices to a file
as binary data or formatted text.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_fwrite</B> <I>(FILE * <VAR>stream</VAR>, const gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX895"></A>
This function writes the elements of the matrix <VAR>m</VAR> to the stream
<VAR>stream</VAR> in binary format.  The return value is 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem writing to the file.  Since the
data is written in the native binary format it may not be portable
between different architectures.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_fread</B> <I>(FILE * <VAR>stream</VAR>, gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX896"></A>
This function reads into the matrix <VAR>m</VAR> from the open stream
<VAR>stream</VAR> in binary format.  The matrix <VAR>m</VAR> must be preallocated
with the correct dimensions since the function uses the size of <VAR>m</VAR> to
determine how many bytes to read.  The return value is 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem reading from the file.  The
data is assumed to have been written in the native binary format on the
same architecture.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_fprintf</B> <I>(FILE * <VAR>stream</VAR>, const gsl_matrix * <VAR>m</VAR>, const char * <VAR>format</VAR>)</I>
<DD><A NAME="IDX897"></A>
This function writes the elements of the matrix <VAR>m</VAR> line-by-line to
the stream <VAR>stream</VAR> using the format specifier <VAR>format</VAR>, which
should be one of the <CODE>%g</CODE>, <CODE>%e</CODE> or <CODE>%f</CODE> formats for
floating point numbers and <CODE>%d</CODE> for integers.  The function returns
0 for success and <CODE>GSL_EFAILED</CODE> if there was a problem writing to
the file.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_fscanf</B> <I>(FILE * <VAR>stream</VAR>, gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX898"></A>
This function reads formatted data from the stream <VAR>stream</VAR> into the
matrix <VAR>m</VAR>.  The matrix <VAR>m</VAR> must be preallocated with the correct
dimensions since the function uses the size of <VAR>m</VAR> to determine how many
numbers to read.  The function returns 0 for success and
<CODE>GSL_EFAILED</CODE> if there was a problem reading from the file.
</DL>

</P>


<H3><A NAME="SEC178" HREF="gsl-ref_toc.html#TOC178">Matrix views</A></H3>

<P>
A matrix view is a temporary object, stored on the stack, which can be
used to operate on a subset of matrix elements.  Matrix views can be
defined for both constant and non-constant matrices using separate types
that preserve constness.  A matrix view has the type
<CODE>gsl_matrix_view</CODE> and a constant matrix view has the type
<CODE>gsl_matrix_const_view</CODE>.  In both cases the elements of the view
can by accessed using the <CODE>matrix</CODE> component of the view object.  A
pointer <CODE>gsl_matrix *</CODE> or <CODE>const gsl_matrix *</CODE> can be obtained
by taking the address of the <CODE>matrix</CODE> component with the <CODE>&#38;</CODE>
operator.  In addition to matrix views it is also possible to create
vector views of a matrix, such as row or column views.

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_matrix_view <B>gsl_matrix_submatrix</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>k1</VAR>, size_t <VAR>k2</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX899"></A>
<DT><U>Function:</U> gsl_matrix_const_view <B>gsl_matrix_const_submatrix</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>k1</VAR>, size_t <VAR>k2</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX900"></A>
These functions return a matrix view of a submatrix of the matrix
<VAR>m</VAR>.  The upper-left element of the submatrix is the element
(<VAR>k1</VAR>,<VAR>k2</VAR>) of the original matrix.  The submatrix has <VAR>n1</VAR>
rows and <VAR>n2</VAR> columns.  The physical number of columns in memory
given by <VAR>tda</VAR> is unchanged.  Mathematically, the
(i,j)-th element of the new matrix is given by,

</P>

<PRE>
m'(i,j) = m-&#62;data[(k1*m-&#62;tda + k2) + i*m-&#62;tda + j]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n1-1</CODE> and the index <VAR>j</VAR>
runs from 0 to <CODE>n2-1</CODE>.

</P>
<P>
The <CODE>data</CODE> pointer of the returned matrix struct is set to null if
the combined parameters (<VAR>i</VAR>,<VAR>j</VAR>,<VAR>n1</VAR>,<VAR>n2</VAR>,<VAR>tda</VAR>)
overrun the ends of the original matrix.

</P>
<P>
The new matrix view is only a view of the block underlying the existing
matrix, <VAR>m</VAR>.  The block containing the elements of <VAR>m</VAR> is not
owned by the new matrix view.  When the view goes out of scope the
original matrix <VAR>m</VAR> and its block will continue to exist.  The
original memory can only be deallocated by freeing the original matrix.
Of course, the original matrix should not be deallocated while the view
is still in use.

</P>
<P>
The function <CODE>gsl_matrix_const_submatrix</CODE> is equivalent to
<CODE>gsl_matrix_submatrix</CODE> but can be used for matrices which are
declared <CODE>const</CODE>.
</DL>

</P>

<P>
<DL>
<DT><U>Function:</U> gsl_matrix_view <B>gsl_matrix_view_array</B> <I>(double * <VAR>base</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX901"></A>
<DT><U>Function:</U> gsl_matrix_const_view <B>gsl_matrix_const_view_array</B> <I>(const double * <VAR>base</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX902"></A>
These functions return a matrix view of the array <VAR>base</VAR>.  The
matrix has <VAR>n1</VAR> rows and <VAR>n2</VAR> columns.  The physical number of
columns in memory is also given by <VAR>n2</VAR>.  Mathematically, the
(i,j)-th element of the new matrix is given by,

</P>

<PRE>
m'(i,j) = base[i*n2 + j]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n1-1</CODE> and the index <VAR>j</VAR>
runs from 0 to <CODE>n2-1</CODE>.

</P>
<P>
The new matrix is only a view of the array <VAR>base</VAR>.  When the view
goes out of scope the original array <VAR>base</VAR> will continue to exist.
The original memory can only be deallocated by freeing the original
array.  Of course, the original array should not be deallocated while
the view is still in use.

</P>
<P>
The function <CODE>gsl_matrix_const_view_array</CODE> is equivalent to
<CODE>gsl_matrix_view_array</CODE> but can be used for matrices which are
declared <CODE>const</CODE>.
</DL>

</P>

<P>
<DL>
<DT><U>Function:</U> gsl_matrix_view <B>gsl_matrix_view_array_with_tda</B> <I>(double * <VAR>base</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>, size_t <VAR>tda</VAR>)</I>
<DD><A NAME="IDX903"></A>
<DT><U>Function:</U> gsl_matrix_const_view <B>gsl_matrix_const_view_array_with_tda</B> <I>(const double * <VAR>base</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>, size_t <VAR>tda</VAR>)</I>
<DD><A NAME="IDX904"></A>
These functions return a matrix view of the array <VAR>base</VAR> with a
physical number of columns <VAR>tda</VAR> which may differ from the corresponding
dimension of the matrix.  The matrix has <VAR>n1</VAR> rows and <VAR>n2</VAR>
columns, and the physical number of columns in memory is given by
<VAR>tda</VAR>.  Mathematically, the (i,j)-th element of the new
matrix is given by,

</P>

<PRE>
m'(i,j) = base[i*tda + j]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n1-1</CODE> and the index <VAR>j</VAR>
runs from 0 to <CODE>n2-1</CODE>.

</P>
<P>
The new matrix is only a view of the array <VAR>base</VAR>.  When the view
goes out of scope the original array <VAR>base</VAR> will continue to exist.
The original memory can only be deallocated by freeing the original
array.  Of course, the original array should not be deallocated while
the view is still in use.

</P>
<P>
The function <CODE>gsl_matrix_const_view_array_with_tda</CODE> is equivalent
to <CODE>gsl_matrix_view_array_with_tda</CODE> but can be used for matrices
which are declared <CODE>const</CODE>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_matrix_view <B>gsl_matrix_view_vector</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX905"></A>
<DT><U>Function:</U> gsl_matrix_const_view <B>gsl_matrix_const_view_vector</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>)</I>
<DD><A NAME="IDX906"></A>
These functions return a matrix view of the vector <VAR>v</VAR>.  The matrix
has <VAR>n1</VAR> rows and <VAR>n2</VAR> columns. The vector must have unit
stride. The physical number of columns in memory is also given by
<VAR>n2</VAR>.  Mathematically, the (i,j)-th element of the new
matrix is given by,

</P>

<PRE>
m'(i,j) = v-&#62;data[i*n2 + j]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n1-1</CODE> and the index <VAR>j</VAR>
runs from 0 to <CODE>n2-1</CODE>.

</P>
<P>
The new matrix is only a view of the vector <VAR>v</VAR>.  When the view
goes out of scope the original vector <VAR>v</VAR> will continue to exist.
The original memory can only be deallocated by freeing the original
vector.  Of course, the original vector should not be deallocated while
the view is still in use.

</P>
<P>
The function <CODE>gsl_matrix_const_view_vector</CODE> is equivalent to
<CODE>gsl_matrix_view_vector</CODE> but can be used for matrices which are
declared <CODE>const</CODE>.
</DL>

</P>

<P>
<DL>
<DT><U>Function:</U> gsl_matrix_view <B>gsl_matrix_view_vector_with_tda</B> <I>(gsl_vector * <VAR>v</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>, size_t <VAR>tda</VAR>)</I>
<DD><A NAME="IDX907"></A>
<DT><U>Function:</U> gsl_matrix_const_view <B>gsl_matrix_const_view_vector_with_tda</B> <I>(const gsl_vector * <VAR>v</VAR>, size_t <VAR>n1</VAR>, size_t <VAR>n2</VAR>, size_t <VAR>tda</VAR>)</I>
<DD><A NAME="IDX908"></A>
These functions return a matrix view of the vector <VAR>v</VAR> with a
physical number of columns <VAR>tda</VAR> which may differ from the
corresponding matrix dimension.  The vector must have unit stride. The
matrix has <VAR>n1</VAR> rows and <VAR>n2</VAR> columns, and the physical number
of columns in memory is given by <VAR>tda</VAR>.  Mathematically, the
(i,j)-th element of the new matrix is given by,

</P>

<PRE>
m'(i,j) = v-&#62;data[i*tda + j]
</PRE>

<P>
where the index <VAR>i</VAR> runs from 0 to <CODE>n1-1</CODE> and the index <VAR>j</VAR>
runs from 0 to <CODE>n2-1</CODE>.

</P>
<P>
The new matrix is only a view of the vector <VAR>v</VAR>.  When the view
goes out of scope the original vector <VAR>v</VAR> will continue to exist.
The original memory can only be deallocated by freeing the original
vector.  Of course, the original vector should not be deallocated while
the view is still in use.

</P>
<P>
The function <CODE>gsl_matrix_const_view_vector_with_tda</CODE> is equivalent
to <CODE>gsl_matrix_view_vector_with_tda</CODE> but can be used for matrices
which are declared <CODE>const</CODE>.
</DL>

</P>



<H3><A NAME="SEC179" HREF="gsl-ref_toc.html#TOC179">Creating row and column views</A></H3>

<P>
In general there are two ways to access an object, by reference or by
copying.  The functions described in this section create vector views
which allow access to a row or column of a matrix by reference.
Modifying elements of the view is equivalent to modifying the matrix,
since both the vector view and the matrix point to the same memory
block.

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_matrix_row</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX909"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_matrix_const_row</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX910"></A>
These functions return a vector view of the <VAR>i</VAR>-th row of the matrix
<VAR>m</VAR>.  The <CODE>data</CODE> pointer of the new vector is set to null if
<VAR>i</VAR> is out of range.

</P>
<P>
The function <CODE>gsl_vector_const_row</CODE> is equivalent to
<CODE>gsl_matrix_row</CODE> but can be used for matrices which are declared
<CODE>const</CODE>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_matrix_column</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX911"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_matrix_const_column</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX912"></A>
These functions return a vector view of the <VAR>j</VAR>-th column of the
matrix <VAR>m</VAR>.  The <CODE>data</CODE> pointer of the new vector is set to
null if <VAR>j</VAR> is out of range.

</P>
<P>
The function <CODE>gsl_vector_const_column</CODE> is equivalent to
<CODE>gsl_matrix_column</CODE> but can be used for matrices which are declared
<CODE>const</CODE>.
</DL>

</P>
<P>
<A NAME="IDX913"></A>
<A NAME="IDX914"></A>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_matrix_diagonal</B> <I>(gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX915"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_matrix_const_diagonal</B> <I>(const gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX916"></A>
These functions returns a vector view of the diagonal of the matrix
<VAR>m</VAR>. The matrix <VAR>m</VAR> is not required to be square. For a
rectangular matrix the length of the diagonal is the same as the smaller
dimension of the matrix.

</P>
<P>
The function <CODE>gsl_matrix_const_diagonal</CODE> is equivalent to
<CODE>gsl_matrix_diagonal</CODE> but can be used for matrices which are
declared <CODE>const</CODE>.
</DL>

</P>
<P>
<A NAME="IDX917"></A>
<A NAME="IDX918"></A>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_matrix_subdiagonal</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>k</VAR>)</I>
<DD><A NAME="IDX919"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_matrix_const_subdiagonal</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>k</VAR>)</I>
<DD><A NAME="IDX920"></A>
These functions return a vector view of the <VAR>k</VAR>-th subdiagonal of
the matrix <VAR>m</VAR>. The matrix <VAR>m</VAR> is not required to be square.
The diagonal of the matrix corresponds to k = 0.

</P>
<P>
The function <CODE>gsl_matrix_const_subdiagonal</CODE> is equivalent to
<CODE>gsl_matrix_subdiagonal</CODE> but can be used for matrices which are
declared <CODE>const</CODE>.
</DL>

</P>
<P>
<A NAME="IDX921"></A>
<A NAME="IDX922"></A>
<DL>
<DT><U>Function:</U> gsl_vector_view <B>gsl_matrix_superdiagonal</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>k</VAR>)</I>
<DD><A NAME="IDX923"></A>
<DT><U>Function:</U> gsl_vector_const_view <B>gsl_matrix_const_superdiagonal</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t <VAR>k</VAR>)</I>
<DD><A NAME="IDX924"></A>
These functions return a vector view of the <VAR>k</VAR>-th superdiagonal of
the matrix <VAR>m</VAR>. The matrix <VAR>m</VAR> is not required to be square. The
diagonal of the matrix corresponds to k = 0.

</P>
<P>
The function <CODE>gsl_matrix_const_superdiagonal</CODE> is equivalent to
<CODE>gsl_matrix_superdiagonal</CODE> but can be used for matrices which are
declared <CODE>const</CODE>.
</DL>

</P>



<H3><A NAME="SEC180" HREF="gsl-ref_toc.html#TOC180">Copying matrices</A></H3>

<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_memcpy</B> <I>(gsl_matrix * <VAR>dest</VAR>, const gsl_matrix * <VAR>src</VAR>)</I>
<DD><A NAME="IDX925"></A>
This function copies the elements of the matrix <VAR>src</VAR> into the
matrix <VAR>dest</VAR>.  The two matrices must have the same size.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_swap</B> <I>(gsl_matrix * <VAR>m1</VAR>, gsl_matrix * <VAR>m2</VAR>)</I>
<DD><A NAME="IDX926"></A>
This function exchanges the elements of the matrices <VAR>m1</VAR> and
<VAR>m2</VAR> by copying.  The two matrices must have the same size.
</DL>

</P>


<H3><A NAME="SEC181" HREF="gsl-ref_toc.html#TOC181">Copying rows and columns</A></H3>

<P>
The functions described in this section copy a row or column of a matrix
into a vector.  This allows the elements of the vector and the matrix to
be modified independently.  Note that if the matrix and the vector point
to overlapping regions of memory then the result will be undefined.  The
same effect can be achieved with more generality using
<CODE>gsl_vector_memcpy</CODE> with vector views of rows and columns.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_get_row</B> <I>(gsl_vector * <VAR>v</VAR>, const gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>)</I>
<DD><A NAME="IDX927"></A>
This function copies the elements of the <VAR>i</VAR>-th row of the matrix
<VAR>m</VAR> into the vector <VAR>v</VAR>.  The length of the vector must be the
same as the length of the row.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_get_col</B> <I>(gsl_vector * <VAR>v</VAR>, const gsl_matrix * <VAR>m</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX928"></A>
This function copies the elements of the <VAR>j</VAR>-th column of the matrix
<VAR>m</VAR> into the vector <VAR>v</VAR>.  The length of the vector must be the
same as the length of the column.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_set_row</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX929"></A>
This function copies the elements of the vector <VAR>v</VAR> into the
<VAR>i</VAR>-th row of the matrix <VAR>m</VAR>.  The length of the vector must be
the same as the length of the row.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_set_col</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>j</VAR>, const gsl_vector * <VAR>v</VAR>)</I>
<DD><A NAME="IDX930"></A>
This function copies the elements of the vector <VAR>v</VAR> into the
<VAR>j</VAR>-th column of the matrix <VAR>m</VAR>.  The length of the vector must be
the same as the length of the column.
</DL>

</P>


<H3><A NAME="SEC182" HREF="gsl-ref_toc.html#TOC182">Exchanging rows and columns</A></H3>

<P>
The following functions can be used to exchange the rows and columns of
a matrix.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_swap_rows</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX931"></A>
This function exchanges the <VAR>i</VAR>-th and <VAR>j</VAR>-th rows of the matrix
<VAR>m</VAR> in-place.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_swap_columns</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX932"></A>
This function exchanges the <VAR>i</VAR>-th and <VAR>j</VAR>-th columns of the
matrix <VAR>m</VAR> in-place.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_swap_rowcol</B> <I>(gsl_matrix * <VAR>m</VAR>, size_t <VAR>i</VAR>, size_t <VAR>j</VAR>)</I>
<DD><A NAME="IDX933"></A>
This function exchanges the <VAR>i</VAR>-th row and <VAR>j</VAR>-th column of the
matrix <VAR>m</VAR> in-place.  The matrix must be square for this operation to
be possible.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_transpose_memcpy</B> <I>(gsl_matrix * <VAR>dest</VAR>, const gsl_matrix * <VAR>src</VAR>)</I>
<DD><A NAME="IDX934"></A>
This function makes the matrix <VAR>dest</VAR> the transpose of the matrix
<VAR>src</VAR> by copying the elements of <VAR>src</VAR> into <VAR>dest</VAR>.  This
function works for all matrices provided that the dimensions of the matrix
<VAR>dest</VAR> match the transposed dimensions of the matrix <VAR>src</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_transpose</B> <I>(gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX935"></A>
This function replaces the matrix <VAR>m</VAR> by its transpose by copying
the elements of the matrix in-place.  The matrix must be square for this
operation to be possible.
</DL>

</P>


<H3><A NAME="SEC183" HREF="gsl-ref_toc.html#TOC183">Matrix operations</A></H3>

<P>
The following operations are defined for real and complex matrices.

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_add</B> <I>(gsl_matrix * <VAR>a</VAR>, const gsl_matrix * <VAR>b</VAR>)</I>
<DD><A NAME="IDX936"></A>
This function adds the elements of matrix <VAR>b</VAR> to the elements of
matrix <VAR>a</VAR>, a'(i,j) = a(i,j) + b(i,j). The two matrices must have the
same dimensions.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_sub</B> <I>(gsl_matrix * <VAR>a</VAR>, const gsl_matrix * <VAR>b</VAR>)</I>
<DD><A NAME="IDX937"></A>
This function subtracts the elements of matrix <VAR>b</VAR> from the elements of
matrix <VAR>a</VAR>, a'(i,j) = a(i,j) - b(i,j). The two matrices must have the
same dimensions.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_mul_elements</B> <I>(gsl_matrix * <VAR>a</VAR>, const gsl_matrix * <VAR>b</VAR>)</I>
<DD><A NAME="IDX938"></A>
This function multiplies the elements of matrix <VAR>a</VAR> by the elements of
matrix <VAR>b</VAR>, a'(i,j) = a(i,j) * b(i,j). The two matrices must have the
same dimensions.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_div_elements</B> <I>(gsl_matrix * <VAR>a</VAR>, const gsl_matrix * <VAR>b</VAR>)</I>
<DD><A NAME="IDX939"></A>
This function divides the elements of matrix <VAR>a</VAR> by the elements of
matrix <VAR>b</VAR>, a'(i,j) = a(i,j) / b(i,j). The two matrices must have the
same dimensions.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_scale</B> <I>(gsl_matrix * <VAR>a</VAR>, const double <VAR>x</VAR>)</I>
<DD><A NAME="IDX940"></A>
This function multiplies the elements of matrix <VAR>a</VAR> by the constant
factor <VAR>x</VAR>, a'(i,j) = x a(i,j).
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_add_constant</B> <I>(gsl_matrix * <VAR>a</VAR>, const double <VAR>x</VAR>)</I>
<DD><A NAME="IDX941"></A>
This function adds the constant value <VAR>x</VAR> to the elements of the
matrix <VAR>a</VAR>, a'(i,j) = a(i,j) + x.
</DL>

</P>


<H3><A NAME="SEC184" HREF="gsl-ref_toc.html#TOC184">Finding maximum and minimum elements of matrices</A></H3>

<P>
The following operations are only defined for real matrices.

</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_matrix_max</B> <I>(const gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX942"></A>
This function returns the maximum value in the matrix <VAR>m</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> double <B>gsl_matrix_min</B> <I>(const gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX943"></A>
This function returns the minimum value in the matrix <VAR>m</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_minmax</B> <I>(const gsl_matrix * <VAR>m</VAR>, double * <VAR>min_out</VAR>, double * <VAR>max_out</VAR>)</I>
<DD><A NAME="IDX944"></A>
This function returns the minimum and maximum values in the matrix
<VAR>m</VAR>, storing them in <VAR>min_out</VAR> and <VAR>max_out</VAR>.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_max_index</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t * <VAR>imax</VAR>, size_t * <VAR>jmax</VAR>)</I>
<DD><A NAME="IDX945"></A>
This function returns the indices of the maximum value in the matrix
<VAR>m</VAR>, storing them in <VAR>imax</VAR> and <VAR>jmax</VAR>.  When there are
several equal maximum elements then the first element found is returned,
searching in row-major order.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_min_index</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t * <VAR>imax</VAR>, size_t * <VAR>jmax</VAR>)</I>
<DD><A NAME="IDX946"></A>
This function returns the indices of the minimum value in the matrix
<VAR>m</VAR>, storing them in <VAR>imax</VAR> and <VAR>jmax</VAR>.  When there are
several equal minimum elements then the first element found is returned,
searching in row-major order.
</DL>

</P>
<P>
<DL>
<DT><U>Function:</U> void <B>gsl_matrix_minmax_index</B> <I>(const gsl_matrix * <VAR>m</VAR>, size_t * <VAR>imin</VAR>, size_t * <VAR>imax</VAR>)</I>
<DD><A NAME="IDX947"></A>
This function returns the indices of the minimum and maximum values in
the matrix <VAR>m</VAR>, storing them in (<VAR>imin</VAR>,<VAR>jmin</VAR>) and
(<VAR>imax</VAR>,<VAR>jmax</VAR>). When there are several equal minimum or maximum
elements then the first elements found are returned, searching in
row-major order.
</DL>

</P>


<H3><A NAME="SEC185" HREF="gsl-ref_toc.html#TOC185">Matrix properties</A></H3>

<P>
<DL>
<DT><U>Function:</U> int <B>gsl_matrix_isnull</B> <I>(const gsl_matrix * <VAR>m</VAR>)</I>
<DD><A NAME="IDX948"></A>
This function returns 1 if all the elements of the matrix <VAR>m</VAR> are
zero, and 0 otherwise.
</DL>

</P>


<H3><A NAME="SEC186" HREF="gsl-ref_toc.html#TOC186">Example programs for matrices</A></H3>

<P>
The program below shows how to allocate, initialize and read from a matrix
using the functions <CODE>gsl_matrix_alloc</CODE>, <CODE>gsl_matrix_set</CODE> and
<CODE>gsl_matrix_get</CODE>.

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_matrix.h&#62;

int
main (void)
{
  int i, j; 
  gsl_matrix * m = gsl_matrix_alloc (10, 3);
  
  for (i = 0; i &#60; 10; i++)
    for (j = 0; j &#60; 3; j++)
      gsl_matrix_set (m, i, j, 0.23 + 100*i + j);
  
  for (i = 0; i &#60; 100; i++)
    for (j = 0; j &#60; 3; j++)
      printf ("m(%d,%d) = %g\n", i, j, 
              gsl_matrix_get (m, i, j));

  return 0;
}
</PRE>

<P>
Here is the output from the program.  The final loop attempts to read
outside the range of the matrix <CODE>m</CODE>, and the error is trapped by
the range-checking code in <CODE>gsl_matrix_get</CODE>.

</P>

<PRE>
$ ./a.out
m(0,0) = 0.23
m(0,1) = 1.23
m(0,2) = 2.23
m(1,0) = 100.23
m(1,1) = 101.23
m(1,2) = 102.23
...
m(9,2) = 902.23
gsl: matrix_source.c:13: ERROR: first index out of range
Default GSL error handler invoked.
Aborted (core dumped)
</PRE>

<P>
The next program shows how to write a matrix to a file.

</P>

<PRE>
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_matrix.h&#62;

int
main (void)
{
  int i, j, k = 0; 
  gsl_matrix * m = gsl_matrix_alloc (100, 100);
  gsl_matrix * a = gsl_matrix_alloc (100, 100);
  
  for (i = 0; i &#60; 100; i++)
    for (j = 0; j &#60; 100; j++)
      gsl_matrix_set (m, i, j, 0.23 + i + j);

  {  
     FILE * f = fopen ("test.dat", "wb");
     gsl_matrix_fwrite (f, m);
     fclose (f);
  }

  {  
     FILE * f = fopen ("test.dat", "rb");
     gsl_matrix_fread (f, a);
     fclose (f);
  }

  for (i = 0; i &#60; 100; i++)
    for (j = 0; j &#60; 100; j++)
      {
        double mij = gsl_matrix_get (m, i, j);
        double aij = gsl_matrix_get (a, i, j);
        if (mij != aij) k++;
      }

  printf ("differences = %d (should be zero)\n", k);
  return (k &#62; 0);
}
</PRE>

<P>
After running this program the file <TT>`test.dat'</TT> should contain the
elements of <CODE>m</CODE>, written in binary format.  The matrix which is read
back in using the function <CODE>gsl_matrix_fread</CODE> should be exactly
equal to the original matrix.

</P>
<P>
The following program demonstrates the use of vector views.  The program
computes the column norms of a matrix.

</P>

<PRE>
#include &#60;math.h&#62;
#include &#60;stdio.h&#62;
#include &#60;gsl/gsl_matrix.h&#62;
#include &#60;gsl/gsl_blas.h&#62;

int
main (void)
{
  size_t i,j;

  gsl_matrix *m = gsl_matrix_alloc (10, 10);

  for (i = 0; i &#60; 10; i++)
    for (j = 0; j &#60; 10; j++)
      gsl_matrix_set (m, i, j, sin (i) + cos (j));

  for (j = 0; j &#60; 10; j++)
    {
      gsl_vector_view column = gsl_matrix_column (m, j);
      double d;

      d = gsl_blas_dnrm2 (&#38;column.vector);

      printf ("matrix column %d, norm = %g\n", j, d);
    }

  gsl_matrix_free (m);

  return 0;
}
</PRE>

<P>
Here is the output of the program, 

</P>

<PRE>
$ ./a.out
matrix column 0, norm = 4.31461
matrix column 1, norm = 3.1205
matrix column 2, norm = 2.19316
matrix column 3, norm = 3.26114
matrix column 4, norm = 2.53416
matrix column 5, norm = 2.57281
matrix column 6, norm = 4.20469
matrix column 7, norm = 3.65202
matrix column 8, norm = 2.08524
matrix column 9, norm = 3.07313
</PRE>

<P>
The results can be confirmed using GNU OCTAVE,

</P>

<PRE>
$ octave
GNU Octave, version 2.0.16.92
octave&#62; m = sin(0:9)' * ones(1,10) 
               + ones(10,1) * cos(0:9); 
octave&#62; sqrt(sum(m.^2))
ans =
  4.3146  3.1205  2.1932  3.2611  2.5342  2.5728
  4.2047  3.6520  2.0852  3.0731
</PRE>



<H2><A NAME="SEC187" HREF="gsl-ref_toc.html#TOC187">References and Further Reading</A></H2>

<P>
The block, vector and matrix objects in GSL follow the <CODE>valarray</CODE>
model of C++.  A description of this model can be found in the following
reference,

</P>

<UL>
<LI>

B. Stroustrup,
<CITE>The C++ Programming Language</CITE> (3rd Ed), 
Section 22.4 Vector Arithmetic.
Addison-Wesley 1997, ISBN 0-201-88954-4.
</UL>

<P><HR><P>
<p>Go to the <A HREF="gsl-ref_1.html">first</A>, <A HREF="gsl-ref_7.html">previous</A>, <A HREF="gsl-ref_9.html">next</A>, <A HREF="gsl-ref_50.html">last</A> section, <A HREF="gsl-ref_toc.html">table of contents</A>.
</BODY>
</HTML>