summaryrefslogtreecommitdiffstats
path: root/ast/chebymap.c
blob: 29f544bf346dd36cb0770955a49c8a21da874569 (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
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
2192
2193
2194
2195
2196
2197
2198
2199
2200
2201
2202
2203
2204
2205
2206
2207
2208
2209
2210
2211
2212
2213
2214
2215
2216
2217
2218
2219
2220
2221
2222
2223
2224
2225
2226
2227
2228
2229
2230
2231
2232
2233
2234
2235
2236
2237
2238
2239
2240
2241
2242
2243
2244
2245
2246
2247
2248
2249
2250
2251
2252
2253
2254
2255
2256
2257
2258
2259
2260
2261
2262
2263
2264
2265
2266
2267
2268
2269
2270
2271
2272
2273
2274
2275
2276
2277
2278
2279
2280
2281
2282
2283
2284
2285
2286
2287
2288
2289
2290
2291
2292
2293
2294
2295
2296
2297
2298
2299
2300
2301
2302
2303
2304
2305
2306
2307
2308
2309
2310
2311
2312
2313
2314
2315
2316
2317
2318
2319
2320
2321
2322
2323
2324
2325
2326
2327
2328
2329
2330
2331
2332
2333
2334
2335
2336
2337
2338
2339
2340
2341
2342
2343
2344
2345
2346
2347
2348
2349
2350
2351
2352
2353
2354
2355
2356
2357
2358
2359
2360
2361
2362
2363
2364
2365
2366
2367
2368
2369
2370
2371
2372
2373
2374
2375
2376
2377
2378
2379
2380
2381
2382
2383
2384
2385
2386
2387
2388
2389
2390
2391
2392
2393
2394
2395
2396
2397
2398
2399
2400
2401
2402
2403
2404
/* To do:

   - what to do about input positions that fall outside the bounding box.
   Have an attribute that can be used to select "set bad" or "extrapolate"?

   - what about overriding astRate ?

   - Providing an iterative inverse requires the Jacobian to be defined.

     for a PolyMap:   if y = C.x^n     then y' = n.C.x^n-1
     for a ChebyMap:  if y = C.Tn(x)   then y' = n.C.Un-1(x)

     where Un-1(x) is the Chebyshev polynomial of the second kind, degree
     (n-1), evaluated at x. Since PolyMap.GetJacobian function uses PolyMaps
     to express the Jacobian of a PolyMap, it would need to use ChebyMaps
     to express the Jacobian of a ChebyMap. But this would mean that
     ChebyMap needs to be able to represent Chebyshev polynomials of the
     second type. In fact the set of powers associated with each coefficient
     would need to indicate somehow whether to use type 1 or type 2 for
     each power. Could use negative powers to indicate type 2, but PolyMap.StoreArrays
     objects to negative powers. StoreArrays could just store them
     without checking, and then call a virtual function to verify the powers
     are OK.

     Simpler for the moment just to disable iterative inverses in
     ChebyMap.

*/


/*
*class++
*  Name:
*     ChebyMap

*  Purpose:
*     Map coordinates using Chebyshev polynomial functions.

*  Constructor Function:
c     astChebyMap
f     AST_CHEBYMAP

*  Description:
*     A ChebyMap is a form of Mapping which performs a Chebyshev polynomial
*     transformation.  Each output coordinate is a linear combination of
*     Chebyshev polynomials of the first kind, of order zero up to a
*     specified maximum order, evaluated at the input coordinates. The
*     coefficients to be used in the linear combination are specified
*     separately for each output coordinate.
*
*     For a 1-dimensional ChebyMap, the forward transformation is defined
*     as follows:
*
*        f(x) = c0.T0(x') + c1.T1(x') + c2.T2(x') + ...
*
*     where:
*        - Tn(x') is the nth Chebyshev polynomial of the first kind:
*             - T0(x') = 1
*             - T1(x') = x'
*             - Tn+1(x') = 2.x'.Tn(x') + Tn-1(x')
*        - x' is the inpux axis value, x, offset and scaled to the range
*          [-1, 1] as x ranges over a specified bounding box, given when the
*          ChebyMap is created. The input positions, x,  supplied to the
*          forward transformation must fall within the bounding box - bad
*          axis values (AST__BAD) are generated for points outside the
*          bounding box.
*
*     For an N-dimensional ChebyMap, the forward transformation is a
*     generalisation of the above form. Each output axis value is the sum
c     of "ncoeff"
f     of NCOEFF
*     terms, where each term is the product of a single coefficient
*     value and N factors of the form Tn(x'_i), where "x'_i" is the
*     normalised value of the i'th input axis value.
*
*     The forward and inverse transformations are defined independantly
*     by separate sets of coefficients, supplied when the ChebyMap is
*     created. If no coefficients are supplied to define the inverse
*     transformation, the
c     astPolyTran
f     AST_POLYTRAN
*     method of the parent PolyMap class can instead be used to create an
*     inverse transformation. The inverse transformation so generated
*     will be a Chebyshev polynomial with coefficients chosen to minimise
*     the residuals left by a round trip (forward transformation followed
*     by inverse transformation).

*  Inheritance:
*     The ChebyMap class inherits from the PolyMap class.

*  Attributes:
*     The ChebyMap class does not define any new attributes beyond those
*     which are applicable to all PolyMaps.

*  Functions:
c     In addition to those functions applicable to all PolyMap, the
c     following functions may also be applied to all ChebyMaps:
f     In addition to those routines applicable to all PolyMap, the
f     following routines may also be applied to all ChebyMaps:
*
c     - astChebyDomain: Get the bounds of the domain of the ChebyMap
f     - AST_CHEBYDOMAIN: Get the bounds of the domain of the ChebyMap

*  Copyright:
*     Copyright (C) 2017 East Asian Observatory.
*     All Rights Reserved.

*  Licence:
*     This program is free software: you can redistribute it and/or
*     modify it under the terms of the GNU Lesser General Public
*     License as published by the Free Software Foundation, either
*     version 3 of the License, or (at your option) any later
*     version.
*
*     This program is distributed in the hope that it will be useful,
*     but WITHOUT ANY WARRANTY; without even the implied warranty of
*     MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
*     GNU Lesser General Public License for more details.
*
*     You should have received a copy of the GNU Lesser General
*     License along with this program.  If not, see
*     <http://www.gnu.org/licenses/>.

*  Authors:
*     DSB: D.S. Berry (EAO)

*  History:
*     1-MAR-2017 (DSB):
*        Original version.
*     30-MAR-2017 (DSB):
*        Over-ride the astFitPoly1DInit and astFitPoly2DInit virtual
*        functions inherited form the PolyMap class.
*class--
*/

/* Module Macros. */
/* ============== */
/* Set the name of the class we are implementing. This indicates to
   the header files that define class interfaces that they should make
   "protected" symbols available. */
#define astCLASS ChebyMap

/* Include files. */
/* ============== */
/* Interface definitions. */
/* ---------------------- */

#include "globals.h"             /* Thread-safe global data access */
#include "error.h"               /* Error reporting facilities */
#include "memory.h"              /* Memory allocation facilities */
#include "object.h"              /* Base Object class */
#include "pointset.h"            /* Sets of points/coordinates */
#include "polymap.h"             /* Polynomial mappings (parent class) */
#include "cmpmap.h"              /* Compound mappings */
#include "chebymap.h"            /* Interface definition for this class */
#include "unitmap.h"             /* Unit mappings */

/* Error code definitions. */
/* ----------------------- */
#include "ast_err.h"             /* AST error codes */

/* C header files. */
/* --------------- */
#include <ctype.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <string.h>
#include <limits.h>
#include <float.h>

/* Module Variables. */
/* ================= */

/* Address of this static variable is used as a unique identifier for
   member of this class. */
static int class_check;

/* Pointers to parent class methods which are extended by this class. */
static int (* parent_getobjsize)( AstObject *, int * );
static int (* parent_equal)( AstObject *, AstObject *, int * );
static void (* parent_polypowers)( AstPolyMap *, double **, int, const int *, double **, int, int, int * );
static AstPolyMap *(*parent_polytran)( AstPolyMap *, int, double, double, int, const double *, const double *, int * );


#ifdef THREAD_SAFE
/* Define how to initialise thread-specific globals. */
#define GLOBAL_inits \
   globals->Class_Init = 0; \

/* Create the function that initialises global data for this module. */
astMAKE_INITGLOBALS(ChebyMap)

/* Define macros for accessing each item of thread specific global data. */
#define class_init astGLOBAL(ChebyMap,Class_Init)
#define class_vtab astGLOBAL(ChebyMap,Class_Vtab)

#include <pthread.h>


#else

/* Define the class virtual function table and its initialisation flag
   as static variables. */
static AstChebyMapVtab class_vtab;   /* Virtual function table */
static int class_init = 0;       /* Virtual function table initialised? */

#endif


/* External Interface Function Prototypes. */
/* ======================================= */
/* The following functions have public prototypes only (i.e. no
   protected prototypes), so we must provide local prototypes for use
   within this module. */
AstChebyMap *astChebyMapId_( int, int, int, const double[], int, const double[],
                             const double[], const double[], const double[],
                             const double[], const char *, ... );


/* Prototypes for Private Member Functions. */
/* ======================================== */
static AstPolyMap *PolyTran( AstPolyMap *, int, double, double, int, const double *, const double *, int * );
static int Equal( AstObject *, AstObject *, int * );
static int GetIterInverse( AstPolyMap *, int * );
static int GetObjSize( AstObject *, int * );
static void ChebyDomain( AstChebyMap *, int, double *, double *, int * );
static void Copy( const AstObject *, AstObject *, int * );
static void Delete( AstObject *obj, int * );
static void Dump( AstObject *, AstChannel *, int * );
static void PolyPowers( AstPolyMap *, double **, int, const int *, double **, int, int, int *);
static void FitPoly1DInit( AstPolyMap *, int, double **, AstMinPackData *, double *, int *);
static void FitPoly2DInit( AstPolyMap *, int, double **, AstMinPackData *, double *, int *);

/* Member functions. */
/* ================= */

static void ChebyDomain( AstChebyMap *this, int forward, double *lbnd,
                         double *ubnd, int *status ){
/*
*++
*  Name:
c     astChebyDomain
f     AST_CHEBYDOMAIN

*  Purpose:
*     Returns the bounding box of the domain of a ChebyMap.

*  Type:
*     Public virtual function.

*  Synopsis:
c     #include "chebymap.h"
c     void astChebyDomain( AstChebyMap *this, int forward, double *lbnd,
c                          double *ubnd )
f     CALL AST_CHEBYDOMAIN( THIS, FORWARD, LBND, UBND, STATUS )

*  Class Membership:
*     ChebyMap method.

*  Description:
c     This function
f     This routine
*     returns the upper and lower limits of the box defining the domain
*     of either the forward or inverse transformation of a ChebyMap. These
*     are the values that were supplied when the ChebyMap was created.

*  Parameters:
c     this
f     THIS = INTEGER (Given)
*        Pointer to the ChebyMap.
c     forward
f     FORWARD = LOGICAL (Given)
c        A non-zero
f        A .TRUE.
*        value indicates that the domain of the ChebyMap's
*        forward transformation is to be returned, while a zero
*        value indicates that the domain of the inverse transformation
*        should be returned.
c     lbnd
f     LBND() = DOUBLE PRECISION (Returned)
c        Pointer to an
f        An
*        array in which to return the lower axis bounds of the ChebyMap
*        domain. The number of elements should be at least equal to the
*        number of ChebyMap inputs (if
c        "forward" is non-zero), or outputs (if "forward" is zero).
f        FORWARD is .TRUE.), or outputs (if FORWARD is .FALSE.).
c     ubnd
f     UBND() = DOUBLE PRECISION (Returned)
c        Pointer to an
f        An
*        array in which to return the upper axis bounds of the ChebyMap
*        domain. The number of elements should be at least equal to the
*        number of ChebyMap inputs (if
c        "forward" is non-zero), or outputs (if "forward" is zero).
f        FORWARD is .TRUE.), or outputs (if FORWARD is .FALSE.).
f     STATUS = INTEGER (Given and Returned)
f        The global status.

*  Notes:
*    - If the requested transformation is undefined (i.e. no
*    transformation coefficients were specified when the ChebyMap was
*    created), this method returns a box determined using the
c    astMapBox
f    AST_MAPBOX
*    method on the opposite transformation, if the opposite
*    transformation is defined.
*    - If the above procedure fails to determine a bounding box, the supplied
*    arrays are filled with AST__BAD values but no error is reported.

*--
*/

/* Local Variables: */
   double *lbnd_o;
   double *offset_o;
   double *offset;
   double *scale_o;
   double *scale;
   double *ubnd_o;
   int fwd_o;
   int iax;
   int nax;
   int nax_o;

/* Check the inherited status. */
   if( !astOK ) return;

/* Get the scales and offsets to use, depending on the value of "forward"
   and whether the ChebyMap has been inverted. */
   if( forward != astGetInvert( this ) ) {
      scale = this->scale_f;
      offset = this->offset_f;
      nax = astGetNin( this );
      scale_o = this->scale_i;
      offset_o = this->offset_i;
      nax_o = astGetNout( this );
      fwd_o = 0;
   } else {
      scale = this->scale_i;
      offset = this->offset_i;
      nax = astGetNout( this );
      scale_o = this->scale_f;
      offset_o = this->offset_f;
      nax_o = astGetNin( this );
      fwd_o = 1;
   }

/* Check the domain is defined. */
   if( scale && offset ) {
      for( iax = 0; iax < nax; iax++ ) {
         if( scale[ iax ] != 0.0 ) {
            lbnd[ iax ] = ( -1.0 - offset[ iax ] ) / scale[ iax ];
            ubnd[ iax ] = ( 1.0 - offset[ iax ] ) / scale[ iax ];
         } else {
            lbnd[ iax ] = AST__BAD;
            ubnd[ iax ] = AST__BAD;
         }
      }

/* If the requested domain is not defined, see if it can be determined
   by transforming the domain of the other transformation into the
   requested input ot putput space. */
   } else if( scale_o && offset_o ){

/* Allocate memory to hold the bounding box in the other space (input or
   output), and then store the bounding box values. */
      lbnd_o = astMalloc( nax_o*sizeof( *lbnd_o ) );
      ubnd_o = astMalloc( nax_o*sizeof( *ubnd_o ) );
      if( astOK ) {
         for( iax = 0; iax < nax_o; iax++ ) {
            if( scale_o[ iax ] != 0.0 ) {
               lbnd_o[ iax ] = ( -1.0 - offset_o[ iax ] ) / scale_o[ iax ];
               ubnd_o[ iax ] = ( 1.0 - offset_o[ iax ] ) / scale_o[ iax ];
            } else {
               lbnd_o[ iax ] = AST__BAD;
               ubnd_o[ iax ] = AST__BAD;
            }
         }

/* Loop round finding the bounds on each input axis of the requested
   transformation. */
         for( iax = 0; iax < nax; iax++ ) {
            astMapBox( this, lbnd_o, ubnd_o, fwd_o, iax, lbnd + iax,
                       ubnd + iax, NULL, NULL );
         }

/* Free resources */
         lbnd_o = astFree( lbnd_o );
         ubnd_o = astFree( ubnd_o );
      }


/* If the domain of the other transformation is not defined, return bad values. */
   } else {
      for( iax = 0; iax < nax; iax++ ) {
         lbnd[ iax ] = AST__BAD;
         ubnd[ iax ] = AST__BAD;
      }
   }
}

static int Equal( AstObject *this_object, AstObject *that_object, int *status ) {
/*
*  Name:
*     Equal

*  Purpose:
*     Test if two ChebyMaps are equivalent.

*  Type:
*     Private function.

*  Synopsis:
*     #include "chebymap.h"
*     int Equal( AstObject *this, AstObject *that, int *status )

*  Class Membership:
*     ChebyMap member function (over-rides the astEqual protected
*     method inherited from the astPolyMap class).

*  Description:
*     This function returns a boolean result (0 or 1) to indicate whether
*     two ChebyMaps are equivalent.

*  Parameters:
*     this
*        Pointer to the first Object (a ChebyMap).
*     that
*        Pointer to the second Object.
*     status
*        Pointer to the inherited status variable.

*  Returned Value:
*     One if the ChebyMaps are equivalent, zero otherwise.

*  Notes:
*     - A value of zero will be returned if this function is invoked
*     with the global status set, or if it should fail for any reason.
*/

/* Local Variables: */
   AstChebyMap *that;
   AstChebyMap *this;
   int i;
   int nin;
   int nout;
   int result;

/* Initialise. */
   result = 0;

/* Check the global error status. */
   if ( !astOK ) return result;

/* Invoke the Equal method inherited from the parent PolyMap class. This
   checks that the PolyMaps are equal. */
   result = (*parent_equal)( this_object, that_object, status );
   if( result ) {

/* Obtain pointers to the two ChebyMap structures. */
      this = (AstChebyMap *) this_object;
      that = (AstChebyMap *) that_object;

/* Check the second object is a ChebyMap. We know the first is a
   ChebyMap since we have arrived at this implementation of the virtual
   function. */
      if( astIsAChebyMap( that ) ) {

/* Get the number of axes in the input bounding box (the original input space). */
         nin = astGetInvert( this ) ? astGetNout( this ) : astGetNin( this );

/* Check the bounding box is the same for both ChebyMaps. */
         if( this->scale_f && that->scale_f &&
             this->offset_f && that->offset_f ) {
            for( i = 0; i < nin && result; i++ ) {
               if( !astEQUAL( this->scale_f[ i ], that->scale_f[ i ] ) ||
                   !astEQUAL( this->offset_f[ i ], that->offset_f[ i ] )){
                  result = 0;
               }
            }
         } else if( this->scale_f || that->scale_f ||
                    this->offset_f || that->offset_f ) {
            result = 0;
         }

/* Get the number of axes in the output bounding box (the original output space). */
         nout = astGetInvert( this ) ? astGetNin( this ) : astGetNout( this );

/* Check the bounding box is the same for both ChebyMaps. */
         if( this->scale_i && that->scale_i &&
             this->offset_i && that->offset_i ) {
            for( i = 0; i < nout && result; i++ ) {
               if( !astEQUAL( this->scale_i[ i ], that->scale_i[ i ] ) ||
                   !astEQUAL( this->offset_i[ i ], that->offset_i[ i ] )){
                  result = 0;
               }
            }
         } else if( this->scale_i || that->scale_i ||
                    this->offset_i || that->offset_i ) {
            result = 0;
         }
      }
   }

/* If an error occurred, clear the result value. */
   if ( !astOK ) result = 0;

/* Return the result, */
   return result;
}

static void FitPoly1DInit( AstPolyMap *this_polymap, int forward, double **table,
                           AstMinPackData *data, double *scales, int *status ){
/*
*  Name:
*     FitPoly1DInit

*  Purpose:
*     Perform initialisation needed for FitPoly1D

*  Type:
*     Private function.

*  Synopsis:
*     #include "chebymap.h"
*     void FitPoly1DInit( AstPolyMap *this, int forward, double **table,
*                         AstMinPackData *data, double *scales, int *status )

*  Class Membership:
*     ChebyMap member function (over-rides the astFitPoly1DInit protected
*     method inherited from the PolyMap class).

*  Description:
*     This function performs initialisation needed for FitPoly1D in the
*     PolyMap class.

*  Parameters:
*     this
*        Pointer to the PolyMap.
*     forward
*        Non-zero if the forward transformation of "this" is being
*        replaced. Zero if the inverse transformation of "this" is being
*        replaced.
*     table
*        Pointer to an array of 2 pointers. Each of these pointers points
*        to an array of "nsamp" doubles, being the scaled and sampled values
*        for x1 and y1 in that order.
*     data
*        Pointer to a structure holding information to pass the the
*        service function invoked by the minimisation function.
*     scales
*        Array holding the scaling factors for the two columns of the table.
*        Multiplying the table values by the scale factor produces PolyMap
*        input or output axis values. The scales are modified on exit to
*        take account of the scaling performed by the ChebyMap Transform
*        method.
*/

/* Local Variables; */
   AstChebyMap *this;
   double *px1;
   double *pxp1;
   double maxx;
   double minx;
   double off;
   double pmax;
   double pmin;
   double scl;
   double x1;
   int k;
   int w1;

/* Check the local error status. */
   if ( !astOK ) return;

/* Get a pointer to the ChebyMap structure. */
   this = (AstChebyMap *) this_polymap;

/* Find the bounds of the supplied x1 values. */
   px1 = table[ 0 ];
   minx = *px1;
   maxx = *px1;
   px1++;
   for( k = 1; k < data->nsamp; k++,px1++ ) {
      if( *px1 > maxx ) {
         maxx = *px1;
      } else if( *px1 < minx ) {
         minx = *px1;
      }
   }

/* Transform the above limits from table values into PolyMap axis values. */
   pmax = maxx*scales[ 0 ];
   pmin = minx*scales[ 0 ];

/* Calculate the scale and offset that map the above bounds onto the range
   [-1,+1], and store them in the ChebyMap. */
   if( pmax != pmin ) {
      scl = 2.0/( pmax - pmin );
      off = -( pmax + pmin )/( pmax - pmin );
   } else if( astOK ){
      astError( AST__BADBX, "astPolyTran(%s): New bounding box has zero width "
                "on axis 1.", status, astGetClass(this));
   }

   if( forward ) {
      this->scale_f = (double *) astFree( this->scale_f );
      this->offset_f = (double *) astFree( this->offset_f );

      this->scale_f = (double *) astMalloc( sizeof( double ) );
      this->offset_f = (double *) astMalloc( sizeof( double ) );
      if( astOK ) {
         this->scale_f[ 0 ] = scl;
         this->offset_f[ 0 ] = off;
      }
   } else {
      this->scale_i = (double *) astFree( this->scale_i );
      this->offset_i = (double *) astFree( this->offset_i );

      this->scale_i = (double *) astMalloc( sizeof( double ) );
      this->offset_i = (double *) astMalloc( sizeof( double ) );
      if( astOK ) {
         this->scale_i[ 0 ] = scl;
         this->offset_i[ 0 ] = off;
      }
   }

/* Get pointers to the supplied x1 values. */
   px1 = table[ 0 ];

/* Get pointers to the location for the next "power" of x1. Here "X to
   the power N" is a metaphor for Tn(x). */
   pxp1 = data->xp1;

/* Loop round all samples. */
   for( k = 0; k < data->nsamp; k++ ) {

/* Get the current x1 value, and scale it into the range [-1,+1]. */
      x1 = *(px1++)*scl*scales[0] + off;

/* Find all the required "powers" of x1 and store them in the "xp1"
   component of the data structure. */
      *(pxp1++) = 1.0;
      *(pxp1++) = x1;
      for( w1 = 2; w1 < data->order; w1++,pxp1++ ) {
         pxp1[ 0 ] = 2.0*x1*pxp1[ -1 ] - pxp1[ -2 ];
      }
   }

/* The scaling representing by the scales[0] value will be performed by
   the astTransform method of the ChebyMap class, so reset teh scales[0]
   value to unity, to avoid the scaling being applied twice. */
   scales[ 0 ] = 1.0;

}

static void FitPoly2DInit( AstPolyMap *this_polymap, int forward, double **table,
                           AstMinPackData *data, double *scales, int *status ){
/*
*  Name:
*     FitPoly2DInit

*  Purpose:
*     Perform initialisation needed for FitPoly2D

*  Type:
*     Private function.

*  Synopsis:
*     #include "chebymap.h"
*     void FitPoly2DInit( AstPolyMap *this, int forward, double **table,
*                         AstMinPackData *data, double *scales, int *status )

*  Class Membership:
*     ChebyMap member function (over-rides the astFitPoly2DInit protected
*     method inherited from the PolyMap class).

*  Description:
*     This function performs initialisation needed for FitPoly2D in the
*     PolyMap class..

*  Parameters:
*     this
*        Pointer to the PolyMap.
*     forward
*        Non-zero if the forward transformation of "this" is being
*        replaced. Zero if the inverse transformation of "this" is being
*        replaced.
*     table
*        Pointer to an array of 4 pointers. Each of these pointers points
*        to an array of "nsamp" doubles, being the scaled and sampled values
*        for x1, x2, y1 or y2 in that order.
*     data
*        Pointer to a structure holding information to pass the the
*        service function invoked by the minimisation function.
*     scales
*        Array holding the scaling factors for the four columns of the table.
*        Multiplying the table values by the scale factor produces PolyMap
*        input or output axis values.
*/

/* Local Variables; */
   AstChebyMap *this;
   double *px1;
   double *px2;
   double *pxp1;
   double *pxp2;
   double maxx;
   double maxy;
   double minx;
   double miny;
   double off[ 2 ];
   double pxmax;
   double pxmin;
   double pymax;
   double pymin;
   double scl[ 2 ];
   double x1;
   double x2;
   int k;
   int w1;
   int w2;

/* Check the local error status. */
   if ( !astOK ) return;

/* Get a pointer to the ChebyMap structure. */
   this = (AstChebyMap *) this_polymap;

/* Find the bounds of the supplied x1 and x2 values. */
   px1 = table[ 0 ];
   px2 = table[ 1 ];
   minx = *px1;
   maxx = *px1;
   miny = *px2;
   maxy = *px2;
   px1++;
   px2++;
   for( k = 1; k < data->nsamp; k++,px1++,px2++ ) {
      if( *px1 > maxx ) {
         maxx = *px1;
      } else if( *px1 < minx ) {
         minx = *px1;
      }
      if( *px2 > maxy ) {
         maxy = *px2;
      } else if( *px2 < miny ) {
         miny = *px2;
      }
   }

/* Transform the above limits from table values into PolyMap axis values. */
   pxmax = maxx*scales[ 0 ];
   pxmin = minx*scales[ 0 ];
   pymax = maxy*scales[ 1 ];
   pymin = miny*scales[ 1 ];

/* Calculate the scale and offset that map the above bounds onto the range
   [-1,+1], and store them in the ChebyMap. */
   if( pxmax != pxmin && pymax != pymin ) {
      scl[ 0 ] = 2.0/( pxmax - pxmin );
      off[ 0 ] = -( pxmax + pxmin )/( pxmax - pxmin );
      scl[ 1 ] = 2.0/( pymax - pymin );
      off[ 1 ] = -( pymax + pymin )/( pymax - pymin );
   } else if( astOK ){
      astError( AST__BADBX, "astPolyTran(%s): New bounding box has zero width "
                "on or both axes.", status, astGetClass(this));
   }

   if( forward ) {
      this->scale_f = (double *) astFree( this->scale_f );
      this->offset_f = (double *) astFree( this->offset_f );

      this->scale_f = (double *) astMalloc( 2*sizeof( double ) );
      this->offset_f = (double *) astMalloc( 2*sizeof( double ) );
      if( astOK ) {
         this->scale_f[ 0 ] = scl[ 0 ];
         this->offset_f[ 0 ] = off[ 0 ];
         this->scale_f[ 1 ] = scl[ 1 ];
         this->offset_f[ 1 ] = off[ 1 ];
      }
   } else {
      this->scale_i = (double *) astFree( this->scale_i );
      this->offset_i = (double *) astFree( this->offset_i );

      this->scale_i = (double *) astMalloc( 2*sizeof( double ) );
      this->offset_i = (double *) astMalloc( 2*sizeof( double ) );
      if( astOK ) {
         this->scale_i[ 0 ] = scl[ 0 ];
         this->offset_i[ 0 ] = off[ 0 ];
         this->scale_i[ 1 ] = scl[ 1 ];
         this->offset_i[ 1 ] = off[ 1 ];
      }
   }

/* Get pointers to the supplied x1 and x2 values. */
   px1 = table[ 0 ];
   px2 = table[ 1 ];

/* Get pointers to the location for the next "power" of x1 anmd x2. Here "X to
   the power N" is a metaphor for Tn(x). */
   pxp1 = data->xp1;
   pxp2 = data->xp2;

/* Loop round all samples. */
   for( k = 0; k < data->nsamp; k++ ) {

/* Get the current x1 and x2 values, and scale them into the range [-1,+1]. */
      x1 = *(px1++)*scl[0]*scales[0] + off[0];
      x2 = *(px2++)*scl[1]*scales[1] + off[1];

/* Find all the required "powers" of x1 and store them in the "xp1"
   component of the data structure. */
      *(pxp1++) = 1.0;
      *(pxp1++) = x1;
      for( w1 = 2; w1 < data->order; w1++,pxp1++ ) {
         pxp1[ 0 ] = 2.0*x1*pxp1[ -1 ] - pxp1[ -2 ];
      }

/* Find all the required "powers" of x2 and store them in the "xp2"
   component of the data structure. */
      *(pxp2++) = 1.0;
      *(pxp2++) = x2;
      for( w2 = 2; w2 < data->order; w2++,pxp2++ ) {
         pxp2[ 0 ] = 2.0*x2*pxp2[ -1 ] - pxp2[ -2 ];
      }
   }

/* The scaling representing by the scales[0] and scales[1] values will be
   performed by the astTransform method of the ChebyMap class, so reset the
   scales[0] and scales[1] values to unity, to avoid the scaling being
   applied twice. */
   scales[ 0 ] = 1.0;
   scales[ 1 ] = 1.0;

}

static int GetIterInverse( AstPolyMap *this, int *status ) {
/*
*  Name:
*     GetIterInverse

*  Purpose:
*     Return the value of the IterInverse attribute.

*  Type:
*     Private function.

*  Synopsis:
*     #include "polymap.h"
*     int GetIterInverse( AstObject *this, int *status )

*  Class Membership:
*     ChebyMap member function (over-rides the astGetIterInverse protected
*     method inherited from the parent PolyMap class).

*  Description:
*     This function returns the value of the IterInverse attribute, which
*     is always zero for a ChebyMap.

*  Parameters:
*     this
*        Pointer to the ChebyMap.
*     status
*        Pointer to the inherited status variable.

*  Returned Value:
*     The IterInverse value.

*  Notes:
*     - A value of zero will be returned if this function is invoked
*     with the global status set, or if it should fail for any reason.
*/

   return 0;
}

static int GetObjSize( AstObject *this_object, int *status ) {
/*
*  Name:
*     GetObjSize

*  Purpose:
*     Return the in-memory size of an Object.

*  Type:
*     Private function.

*  Synopsis:
*     #include "chebymap.h"
*     int GetObjSize( AstObject *this, int *status )

*  Class Membership:
*     ChebyMap member function (over-rides the astGetObjSize protected
*     method inherited from the parent class).

*  Description:
*     This function returns the in-memory size of the supplied ChebyMap,
*     in bytes.

*  Parameters:
*     this
*        Pointer to the ChebyMap.
*     status
*        Pointer to the inherited status variable.

*  Returned Value:
*     The Object size, in bytes.

*  Notes:
*     - A value of zero will be returned if this function is invoked
*     with the global status set, or if it should fail for any reason.
*/

/* Local Variables: */
   AstChebyMap *this;
   int nin;
   int nout;
   int result;

/* Initialise. */
   result = 0;

/* Check the global error status. */
   if ( !astOK ) return result;

/* Obtain a pointers to the ChebyMap structure. */
   this = (AstChebyMap *) this_object;

/* Get the number of input and output axes. */
   nin = astGetInvert( this ) ? astGetNout( this ) : astGetNin( this );
   nout = astGetInvert( this ) ? astGetNin( this ) : astGetNout( this );

/* Invoke the GetObjSize method inherited from the parent class, and then
   add on any components of the class structure defined by this class
   which are stored in dynamically allocated memory. */
   result = (*parent_getobjsize)( this_object, status );

   if( this->scale_f ) result += nin*sizeof( *this->scale_f );
   if( this->offset_f ) result += nin*sizeof( *this->offset_f );
   if( this->scale_i ) result += nout*sizeof( *this->scale_i );
   if( this->offset_i ) result += nout*sizeof( *this->offset_i );

/* If an error occurred, clear the result value. */
   if ( !astOK ) result = 0;

/* Return the result, */
   return result;
}

void astInitChebyMapVtab_(  AstChebyMapVtab *vtab, const char *name, int *status ) {
/*
*+
*  Name:
*     astInitChebyMapVtab

*  Purpose:
*     Initialise a virtual function table for a ChebyMap.

*  Type:
*     Protected function.

*  Synopsis:
*     #include "chebymap.h"
*     void astInitChebyMapVtab( AstChebyMapVtab *vtab, const char *name )

*  Class Membership:
*     ChebyMap vtab initialiser.

*  Description:
*     This function initialises the component of a virtual function
*     table which is used by the ChebyMap class.

*  Parameters:
*     vtab
*        Pointer to the virtual function table. The components used by
*        all ancestral classes will be initialised if they have not already
*        been initialised.
*     name
*        Pointer to a constant null-terminated character string which contains
*        the name of the class to which the virtual function table belongs (it
*        is this pointer value that will subsequently be returned by the Object
*        astClass function).
*-
*/

/* Local Variables: */
   astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
   AstObjectVtab *object;        /* Pointer to Object component of Vtab */
   AstPolyMapVtab *polymap;      /* Pointer to PolyMap component of Vtab */

/* Check the local error status. */
   if ( !astOK ) return;

/* Get a pointer to the thread specific global data structure. */
   astGET_GLOBALS(NULL);

/* Initialize the component of the virtual function table used by the
   parent class. */
   astInitPolyMapVtab( (AstPolyMapVtab *) vtab, name );

/* Store a unique "magic" value in the virtual function table. This
   will be used (by astIsAChebyMap) to determine if an object belongs
   to this class.  We can conveniently use the address of the (static)
   class_check variable to generate this unique value. */
   vtab->id.check = &class_check;
   vtab->id.parent = &(((AstPolyMapVtab *) vtab)->id);

/* Initialise member function pointers. */
/* ------------------------------------ */
/* Store pointers to the member functions (implemented here) that provide
   virtual methods for this class. */
/* none */

/* Save the inherited pointers to methods that will be extended, and
   replace them with pointers to the new member functions. */
   object = (AstObjectVtab *) vtab;
   polymap = (AstPolyMapVtab *) vtab;

   polymap->GetIterInverse = GetIterInverse;

   parent_getobjsize = object->GetObjSize;
   object->GetObjSize = GetObjSize;

   parent_polypowers = polymap->PolyPowers;
   polymap->PolyPowers = PolyPowers;

   parent_polytran = polymap->PolyTran;
   polymap->PolyTran = PolyTran;

   parent_equal = object->Equal;
   object->Equal = Equal;

   polymap->FitPoly1DInit = FitPoly1DInit;
   polymap->FitPoly2DInit = FitPoly2DInit;

/* Store pointers to the member functions (implemented here) that
   provide virtual methods for this class. */
   vtab->ChebyDomain = ChebyDomain;

/* Declare the destructor and copy constructor. */
   astSetDelete( (AstObjectVtab *) vtab, Delete );
   astSetCopy( (AstObjectVtab *) vtab, Copy );

/* Declare the class dump function. */
   astSetDump( vtab, Dump, "ChebyMap", "Chebyshev polynomial transformation" );

/* If we have just initialised the vtab for the current class, indicate
   that the vtab is now initialised, and store a pointer to the class
   identifier in the base "object" level of the vtab. */
   if( vtab == &class_vtab ) {
      class_init = 1;
      astSetVtabClassIdentifier( vtab, &(vtab->id) );
   }
}

static void PolyPowers( AstPolyMap *this_polymap, double **work, int ncoord,
                        const int *mxpow, double **ptr, int point, int fwd,
                        int *status ){
/*
*  Name:
*     PolyPowers

*  Purpose:
*     Find the required powers of the input axis values.

*  Type:
*     Private function.

*  Synopsis:
*     #include "chebymap.h"
*     void PolyPowers( AstPolyMap *this, double **work, int ncoord,
*                      const int *mxpow, double **ptr, int point,
*                      int fwd, int *status )

*  Class Membership:
*     ChebyMap member function (over-rides the astPolyPowers protected
*     method inherited from the PolyMap class).

*  Description:
*     This function is used by astTransform to calculate the powers of
*     the axis values for a single input position. In the case of
*     sub-classes, the powers may not be simply powers of the supplied
*     axis values but may be more complex quantities such as a Chebyshev
*     polynomial of the required degree evaluated at the input axis values.

*  Parameters:
*     this
*        Pointer to the PolyMap.
*     work
*        An array of "ncoord" pointers, each pointing to an array of
*        length "max(2,mxpow)". The required values are placed in this
*        array on exit.
*     ncoord
*        The number of axes.
*     mxpow
*        Pointer to an array holding the maximum power required of each
*        axis value. Should have "ncoord" elements.
*     ptr
*        An array of "ncoord" pointers, each pointing to an array holding
*        the axis values. Each of these arrays of axis values must have
*        at least "point+1" elements.
*     point
*        The zero based index of the point within "ptr" that holds the
*        axis values to be exponentiated.
*     fwd
*        Do the supplied coefficients define the foward transformation of
*        the PolyMap?
*/

/* Local Variables; */
   AstChebyMap *this;
   double *scales;
   double *offsets;
   double *pwork;
   double *t;
   double x;
   int coord;
   int ip;

/* Check the local error status. */
   if ( !astOK ) return;

/* Get a pointer to the ChebyMap structure. */
   this = (AstChebyMap *) this_polymap;

/* Either transformation of a ChebyMap (forward or inverse) can be
   defined either as Chebyshev polynomial or as a standard polynomial.
   Chebyshev polynomials always have non-NULL scale array pointers.
   If the scale array pointer is NULL, then the transformation is a
   standard polynomial. If the coefficients relate to a standard
   polynomial, then invoke the astPolyPowers implementation of the parent
   class (PolyMap). */
   if( (fwd && !this->scale_f) || (!fwd && !this->scale_i) ) {
      (*parent_polypowers)( this_polymap, work, ncoord, mxpow, ptr, point,
                            fwd, status );

/* If the coefficients relate to a Chebyshev polynomial... */
   } else {
      scales = fwd ? this->scale_f : this->scale_i;
      offsets = fwd ? this->offset_f : this->offset_i;

/* This method uses a Chebyshev polynomial of the first kind of degree "i"
   evaluated at "x'" instead of "x raised to the power i". Here, "x'" is
   the input axis value scaled and shifted into the range [-1,+1] on each
   axis. Loop over all input axes. */
      for( coord = 0; coord < ncoord; coord++ ) {

/* Get a pointer to the array in which the powers of the current axis
   value are to be returned. */
         pwork = work[ coord ];

/* The Chebyshev function (type 1) of degree zero is always 1.0, regardless
   of the value of x. */
         pwork[ 0 ] = 1.0;

/* Get the input axis value. If it is bad, store bad values for all
   remaining powers. */
         x = ptr[ coord ][ point ];
         if( x == AST__BAD ) {
            for( ip = 1; ip <= mxpow[ coord ]; ip++ ) pwork[ ip ] = AST__BAD;

/* Otherwise, apply the required scaling to the input */
         } else {
            x = x*scales[ coord ] + offsets[ coord ];

/* Return bad values for input positions outside the bounding box
   associated with the transformation. */
            if( fabs( x ) <= 1.0 ) {

/* The Chebyshev function of degree one is equal to x. */
               t = pwork + 1;
               *t = x;

/* Form and store the remaining Chebyshev polynomial values at the input axis value.
   Use the standard recurrence relation: Tn+1(x') = 2.x'.Tn(x') + Tn-1(x'). */
               for( ip = 2; ip <= mxpow[ coord ]; ip++,t++ ) {
                  t[ 1 ] = 2.0*x*t[ 0 ] - t[ -1 ];
               }
            } else {
               for( ip = 1; ip <= mxpow[ coord ]; ip++ ) pwork[ ip ] = AST__BAD;
            }
         }
      }
   }
}

static AstPolyMap *PolyTran( AstPolyMap *this_polymap, int forward, double acc,
                             double maxacc, int maxorder, const double *lbnd,
                             const double *ubnd, int *status ){
/*
*  Name:
*     PolyTran

*  Purpose:
*     Fit a PolyMap inverse or forward transformation.

*  Type:
*     Private function.

*  Synopsis:
*     #include "polymap.h"
*     AstPolyMap *PolyTran( AstPolyMap *this, int forward, double acc,
*                           double maxacc, int maxorder, const double *lbnd,
*                           const double *ubnd )

*  Class Membership:
*     ChebyMap member function (over-rides the astPolyTran method inherited
*     from the PolyMap class).

*  Description:
*     This function creates a new PolyMap which is a copy of the supplied
*     PolyMap, in which a specified transformation (forward or inverse)
*     has been replaced by a new polynomial transformation. The
*     coefficients of the new transformation are estimated by sampling
*     the other transformation and performing a least squares polynomial
*     fit in the opposite direction to the sampled positions and values.
*
*     This method can only be used on (1-input,1-output) or (2-input,2-output)
*     PolyMaps.
*
*     The transformation to create is specified by the "forward" parameter.
*     In what follows "X" refers to the inputs of the PolyMap, and "Y" to
*     the outputs of the PolyMap. The forward transformation transforms
*     input values (X) into output values (Y), and the inverse transformation
*     transforms output values (Y) into input values (X). Within a PolyMap,
*     each transformation is represented by an independent set of
*     polynomials, P_f or P_i: Y=P_f(X) for the forward transformation and
*     X=P_i(Y) for the inverse transformation.
*
*     The "forward" parameter specifies the transformation to be replaced. If
*     it is non-zero, a new forward transformation is created by first finding
*     the input values (X) using the inverse transformation
*     (which must be available) at a regular grid of points (Y) covering a
*     rectangular region of the PolyMap's output space. The coefficients of
*     the required forward polynomial, Y=P_f(X), are chosen in order to
*     minimise the sum of the squared residuals between the sampled values
*     of Y and P_f(X).
*
*     If "forward" is zero (probably the most likely case),
*     a new inverse transformation is created by
*     first finding the output values (Y) using the forward transformation
*     (which must be available) at a regular grid of points (X) covering a
*     rectangular region of the PolyMap's input space. The coefficients of
*     the required inverse polynomial, X=P_i(Y), are chosen in order to
*     minimise the sum of the squared residuals between the sampled values
*     of X and P_i(Y).
*
*     This fitting process is performed repeatedly with increasing
*     polynomial orders (starting with linear) until the target
*     accuracy is achieved, or a specified maximum order is reached. If
*     the target accuracy cannot be achieved even with this maximum-order
*     polynomial, the best fitting maximum-order polynomial is returned so
*     long as its accuracy is better than "maxacc".
*     If it is not, an error is reported.

*  Parameters:
*     this
*        Pointer to the original Mapping.
*     forward
*        If non-zero, the forward PolyMap transformation is replaced.
*        Otherwise the inverse transformation is replaced.
*     acc
*        The target accuracy, expressed as a geodesic distance within
*        the PolyMap's input space (if "forward" is zero)  or output
*        space (if "forward" is non-zero).
*     maxacc
*        The maximum allowed accuracy for an acceptable polynomial,
*        expressed as a geodesic distance within the PolyMap's input
*        space (if "forward" is zero)  or output space (if "forward" is
*        non-zero).
*     maxorder
*        The maximum allowed polynomial order. This is one more than the
*        maximum power of either input axis. So for instance, a value of
*        3 refers to a quadratic polynomial. Note, cross terms with total
*        powers greater than or equal to maxorder are not inlcuded in the
*        fit. So the maximum number of terms in each of the fitted
*        polynomials is maxorder*(maxorder+1)/2.
*     lbnd
*        Pointer to an array holding the lower bounds of a rectangular
*        region within the PolyMap's input space (if "forward" is zero)
*        or output space (if "forward" is non-zero). The new polynomial
*        will be evaluated over this rectangle. The length of this array
*        should equal the value of the PolyMap's Nin or Nout attribute,
*        depending on "forward". If a NULL pointer is supplied, the lower
*        bounds of the box supplied when the ChebyMap was constructed is
*        used.
*     ubnd
*        Pointer to an
*        array holding the upper bounds of a rectangular region within
*        the PolyMap's input space (if "forward" is zero)  or output space
*        (if "forward" is non-zero). The new polynomial will be evaluated
*        over this rectangle.  The length of this array should equal the
*        value of the PolyMap's Nin or Nout attribute, depending on "forward".
*        If a NULL pointer is supplied, the upper bounds of the box supplied
*        when the ChebyMap was constructed is used.

*  Returned Value:
*     astPolyTran()
*        A pointer to the new PolyMap. A NULL pointer will be returned if
*        the fit fails to achieve the accuracy specified by "maxacc", but
*        no error will be reported.

*  Notes:
*     - This function can only be used on 1D or 2D PolyMaps which have
*     the same number of inputs and outputs.
*     - A null Object pointer (AST__NULL) will be returned if this
*     function is invoked with the AST error status set, or if it
*     should fail for any reason.
*/

/* Local Variables: */
   AstChebyMap *this;
   AstPolyMap *result;
   const char *word;
   double *offset;
   double *scale;
   double this_lbnd[ 2 ];
   double this_ubnd[ 2 ];
   int inverted;
   int nax;

/* Initialise. */
   result = NULL;

/* Check the inherited status. */
   if ( !astOK ) return result;

/* Get a pointer to the CHebyMap structure. */
   this = (AstChebyMap *) this_polymap;

/* Select the ChebyMap scales and offsets to be used. */
   inverted = astGetInvert( this );
   if( ( inverted && !forward ) || ( !inverted && forward ) ) {
      word = "inverse";
      scale = this->scale_i;
      offset = this->offset_i;
      nax = ((AstMapping *)this)->nout;
   } else {
      word = "forward";
      scale = this->scale_f;
      offset = this->offset_f;
      nax = ((AstMapping *)this)->nin;
   }

/* The scaled box for a Chebyshev polynomial spans [-1,+1] on each axis.
   Create the corresponding unscaled box. If the user supplies any bounds,
   use them in preference to the bounds in the ChebyMap. */
   if( lbnd ) {
      this_lbnd[ 0 ] = lbnd[ 0 ];
      if( nax > 1 ) this_lbnd[ 1 ] = lbnd[ 1 ];

   } else if( scale && offset ) {
      this_lbnd[ 0 ] = ( -1.0 - offset[ 0 ] )/scale[ 0 ];
      if( nax > 1 ) this_lbnd[ 1 ] = ( -1.0 - offset[ 1 ] )/scale[ 1 ];

   } else if( astOK ) {
      astError( AST__NOBOX, "astPolyTran(%s): The %s transformation is "
                "not a Chebyshev polynomial and therefore requires a "
                "user-supplied bounding box. But no lower bounds were "
                "supplied. ", status, astGetClass( this ), word );
   }


   if( ubnd ) {
      this_ubnd[ 0 ] = ubnd[ 0 ];
      if( nax > 1 ) this_ubnd[ 1 ] = ubnd[ 1 ];
   } else if( scale && offset ) {
      this_ubnd[ 0 ] = ( 1.0 - offset[ 0 ] )/scale[ 0 ];
      if( nax > 1 ) this_ubnd[ 1 ] = ( 1.0 - offset[ 1 ] )/scale[ 1 ];
   } else if( astOK ) {
      astError( AST__NOBOX, "astPolyTran(%s): The %s transformation is "
                "not a Chebyshev polynomial and therefore requires a "
                "user-supplied bounding box. But no upper bounds were "
                "supplied. ", status, astGetClass( this ), word );
   }

/* Invoke the parent astPolyMap method, using the bounding box selected
   above. */
   result = (*parent_polytran)( this_polymap, forward, acc, maxacc, maxorder,
                                this_lbnd, this_ubnd, status );

/* Return the new ChebyMap. */
   return result;
}


/* Functions which access class attributes. */
/* ---------------------------------------- */
/* Implement member functions to access the attributes associated with
   this class using the macros defined for this purpose in the
   "object.h" file. For a description of each attribute, see the class
   interface (in the associated .h file). */

/* Copy constructor. */
/* ----------------- */
static void Copy( const AstObject *objin, AstObject *objout, int *status ) {
/*
*  Name:
*     Copy

*  Purpose:
*     Copy constructor for ChebyMap objects.

*  Type:
*     Private function.

*  Synopsis:
*     void Copy( const AstObject *objin, AstObject *objout, int *status )

*  Description:
*     This function implements the copy constructor for ChebyMap objects.

*  Parameters:
*     objin
*        Pointer to the object to be copied.
*     objout
*        Pointer to the object being constructed.
*     status
*        Pointer to the inherited status variable.

*  Returned Value:
*     void

*  Notes:
*     -  This constructor makes a deep copy, including a copy of the
*     coefficients associated with the input ChebyMap.
*/


/* Local Variables: */
   AstChebyMap *in;               /* Pointer to input ChebyMap */
   AstChebyMap *out;              /* Pointer to output ChebyMap */
   int nin;                       /* No. of input coordinates */
   int nout;                      /* No. of output coordinates */

/* Check the global error status. */
   if ( !astOK ) return;

/* Obtain pointers to the input and output ChebyMaps. */
   in = (AstChebyMap *) objin;
   out = (AstChebyMap *) objout;

/* Nullify the pointers stored in the output object since these will
   currently be pointing at the input data (since the output is a simple
   byte-for-byte copy of the input). Otherwise, the input data could be
   freed by accidient if the output object is deleted due to an error
   occuring in this function. */
   out->scale_f = NULL;
   out->offset_f = NULL;
   out->scale_i = NULL;
   out->offset_i = NULL;

/* Get the number of inputs and outputs of the uninverted Mapping. */
   nin = ( (AstMapping *) in )->nin;
   nout = ( (AstMapping *) in )->nout;

/* Copy the bounding box arrays. */
   if( in->scale_f ) out->scale_f = (double *) astStore( NULL,
                                       (void *) in->scale_f,
                                       sizeof( double )*nin );
   if( in->offset_f ) out->offset_f = (double *) astStore( NULL,
                                       (void *) in->offset_f,
                                       sizeof( double )*nin );
   if( in->scale_i ) out->scale_i = (double *) astStore( NULL,
                                       (void *) in->scale_i,
                                       sizeof( double )*nout );
   if( in->offset_i ) out->offset_i = (double *) astStore( NULL,
                                       (void *) in->offset_i,
                                       sizeof( double )*nout );
}

/* Destructor. */
/* ----------- */
static void Delete( AstObject *obj, int *status ) {
/*
*  Name:
*     Delete

*  Purpose:
*     Destructor for ChebyMap objects.

*  Type:
*     Private function.

*  Synopsis:
*     void Delete( AstObject *obj, int *status )

*  Description:
*     This function implements the destructor for ChebyMap objects.

*  Parameters:
*     obj
*        Pointer to the object to be deleted.
*     status
*        Pointer to the inherited status variable.

*  Returned Value:
*     void

*  Notes:
*     This function attempts to execute even if the global error status is
*     set.
*/

/* Local Variables: */
   AstChebyMap *this;

/* Obtain a pointer to the ChebyMap structure. */
   this = (AstChebyMap *) obj;

/* Free the boundib box arrays. */
   this->scale_f = astFree( this->scale_f );
   this->offset_f = astFree( this->offset_f );
   this->scale_i = astFree( this->scale_i );
   this->offset_i = astFree( this->offset_i );
}

/* Dump function. */
/* -------------- */
static void Dump( AstObject *this_object, AstChannel *channel, int *status ) {
/*
*  Name:
*     Dump

*  Purpose:
*     Dump function for ChebyMap objects.

*  Type:
*     Private function.

*  Synopsis:
*     void Dump( AstObject *this, AstChannel *channel, int *status )

*  Description:
*     This function implements the Dump function which writes out data
*     for the ChebyMap class to an output Channel.

*  Parameters:
*     this
*        Pointer to the ChebyMap whose data are being written.
*     channel
*        Pointer to the Channel to which the data are being written.
*     status
*        Pointer to the inherited status variable.
*/

#define KEY_LEN 50               /* Maximum length of a keyword */

/* Local Variables: */
   AstChebyMap *this;             /* Pointer to the ChebyMap structure */
   char buff[ KEY_LEN + 1 ];     /* Buffer for keyword string */
   char comm[ 100 ];             /* Buffer for comment string */
   int i;                        /* Loop index */
   int nin;                      /* No. of input coords */
   int nout;                     /* No. of output coords */

/* Check the global error status. */
   if ( !astOK ) return;

/* Obtain a pointer to the ChebyMap structure. */
   this = (AstChebyMap *) this_object;

/* Find the number of inputs and outputs of the uninverted Mapping. */
   nin = ( (AstMapping *) this )->nin;
   nout = ( (AstMapping *) this )->nout;

/* Write out values representing the instance variables for the
   ChebyMap class.  */

/* The input axis scale factors. */
   if( this->scale_f ){
      for( i = 0; i < nin; i++ ){
         (void) sprintf( buff, "FSCL%d", i + 1 );
         (void) sprintf( comm, "Scale factor on input %d", i + 1 );
         astWriteDouble( channel, buff, 1, 1, (this->scale_f)[ i ], comm );
      }
   }

/* The input axis offsets. */
   if( this->offset_f ){
      for( i = 0; i < nin; i++ ){
         (void) sprintf( buff, "FOFF%d", i + 1 );
         (void) sprintf( comm, "Offset on input %d", i + 1 );
         astWriteDouble( channel, buff, 1, 1, (this->offset_f)[ i ], comm );
      }
   }

/* The output axis scale factors. */
   if( this->scale_i ){
      for( i = 0; i < nout; i++ ){
         (void) sprintf( buff, "ISCL%d", i + 1 );
         (void) sprintf( comm, "Scale factor on output %d", i + 1 );
         astWriteDouble( channel, buff, 1, 1, (this->scale_i)[ i ], comm );
      }
   }

/* The output axis offsets. */
   if( this->offset_i ){
      for( i = 0; i < nout; i++ ){
         (void) sprintf( buff, "IOFF%d", i + 1 );
         (void) sprintf( comm, "Offset on output %d", i + 1 );
         astWriteDouble( channel, buff, 1, 1, (this->offset_i)[ i ], comm );
      }
   }

/* Undefine macros local to this function. */
#undef KEY_LEN
}

/* Standard class functions. */
/* ========================= */
/* Implement the astIsAChebyMap and astCheckChebyMap functions using the macros
   defined for this purpose in the "object.h" header file. */
astMAKE_ISA(ChebyMap,Mapping)
astMAKE_CHECK(ChebyMap)

AstChebyMap *astChebyMap_( int nin, int nout, int ncoeff_f, const double coeff_f[],
                           int ncoeff_i, const double coeff_i[],
                           const double lbnd_f[], const double ubnd_f[],
                           const double lbnd_i[], const double ubnd_i[],
                           const char *options, int *status, ...){
/*
*++
*  Name:
c     astChebyMap
f     AST_CHEBYMAP

*  Purpose:
*     Create a ChebyMap.

*  Type:
*     Public function.

*  Synopsis:
c     #include "chebymap.h"
c     AstChebyMap *astChebyMap( int nin, int nout, int ncoeff_f, const double coeff_f[],
c                               int ncoeff_i, const double coeff_i[],
c                               const double lbnd_f[], const double ubnd_f[],
c                               const double lbnd_i[], const double ubnd_i[],
c                               const char *options, ... )
f     RESULT = AST_CHEBYMAP( NIN, NOUT, NCOEFF_F, COEFF_F, NCOEFF_I, COEFF_I,
f                            LBND_F, UBND_F, LBND_I, UBND_I, OPTIONS, STATUS )

*  Class Membership:
*     ChebyMap constructor.

*  Description:
*     This function creates a new ChebyMap and optionally initialises
*     its attributes.
*
*     A ChebyMap is a form of Mapping which performs a Chebyshev polynomial
*     transformation.  Each output coordinate is a linear combination of
*     Chebyshev polynomials of the first kind, of order zero up to a
*     specified maximum order, evaluated at the input coordinates. The
*     coefficients to be used in the linear combination are specified
*     separately for each output coordinate.
*
*     For a 1-dimensional ChebyMap, the forward transformation is defined
*     as follows:
*
*        f(x) = c0.T0(x') + c1.T1(x') + c2.T2(x') + ...
*
*     where:
*        - Tn(x') is the nth Chebyshev polynomial of the first kind:
*             - T0(x') = 1
*             - T1(x') = x'
*             - Tn+1(x') = 2.x'.Tn(x') + Tn-1(x')
*        - x' is the inpux axis value, x, offset and scaled to the range
*          [-1, 1] as x ranges over a specified bounding box, given when the
*          ChebyMap is created. The input positions, x,  supplied to the
*          forward transformation must fall within the bounding box - bad
*          axis values (AST__BAD) are generated for points outside the
*          bounding box.
*
*     For an N-dimensional ChebyMap, the forward transformation is a
*     generalisation of the above form. Each output axis value is the sum
c     of "ncoeff"
f     of NCOEFF
*     terms, where each term is the product of a single coefficient
*     value and N factors of the form Tn(x'_i), where "x'_i" is the
*     normalised value of the i'th input axis value.
*
*     The forward and inverse transformations are defined independantly
*     by separate sets of coefficients, supplied when the ChebyMap is
*     created. If no coefficients are supplied to define the inverse
*     transformation, the
c     astPolyTran
f     AST_POLYTRAN
*     method of the parent PolyMap class can instead be used to create an
*     inverse transformation. The inverse transformation so generated
*     will be a Chebyshev polynomial with coefficients chosen to minimise
*     the residuals left by a round trip (forward transformation followed
*     by inverse transformation).

*  Parameters:
c     nin
f     NIN = INTEGER (Given)
*        The number of input coordinates.
c     nout
f     NOUT = INTEGER (Given)
*        The number of output coordinates.
c     ncoeff_f
f     NCOEFF_F = INTEGER (Given)
*        The number of non-zero coefficients necessary to define the
*        forward transformation of the ChebyMap. If zero is supplied, the
*        forward transformation will be undefined.
c     coeff_f
f     COEFF_F( * ) = DOUBLE PRECISION (Given)
*        An array containing
c        "ncoeff_f*( 2 + nin )" elements. Each group of "2 + nin"
f        "NCOEFF_F*( 2 + NIN )" elements. Each group of "2 + NIN"
*        adjacent elements describe a single coefficient of the forward
*        transformation. Within each such group, the first element is the
*        coefficient value; the next element is the integer index of the
*        ChebyMap output which uses the coefficient within its defining
*        expression (the first output has index 1); the remaining elements
*        of the group give the integer powers to use with each input
*        coordinate value (powers must not be negative, and floating
*        point values are rounded to the nearest integer).
c        If "ncoeff_f" is zero, a NULL pointer may be supplied for "coeff_f".
*
*        For instance, if the ChebyMap has 3 inputs and 2 outputs, each group
*        consisting of 5 elements, A groups such as "(1.2, 2.0, 1.0, 3.0, 0.0)"
*        describes a coefficient with value 1.2 which is used within the
*        definition of output 2. The output value is incremented by the
*        product of the coefficient value, the value of the Chebyshev
*        polynomial of power 1 evaluated at input coordinate 1, and the
*        value of the Chebyshev polynomial of power 3 evaluated at input
*        coordinate 2. Input coordinate 3 is not used since its power is
*        specified as zero. As another example, the group "(-1.0, 1.0,
*        0.0, 0.0, 0.0 )" adds a constant value -1.0 onto output 1 (it is
*        a constant value since the power for every input axis is given as
*        zero).
*
c        Each final output coordinate value is the sum of the "ncoeff_f" terms
c        described by the "ncoeff_f" groups within the supplied array.
f        Each final output coordinate value is the sum of the "NCOEFF_F" terms
f        described by the "NCOEFF_F" groups within the supplied array.
c     ncoeff_i
f     NCOEFF_I = INTEGER (Given)
*        The number of non-zero coefficients necessary to define the
*        inverse transformation of the ChebyMap. If zero is supplied, the
*        inverse transformation will be undefined.
c     coeff_i
f     COEFF_I( * ) = DOUBLE PRECISION (Given)
*        An array containing
c        "ncoeff_i*( 2 + nout )" elements. Each group of "2 + nout"
f        "NCOEFF_I*( 2 + NOUT )" elements. Each group of "2 + NOUT"
*        adjacent elements describe a single coefficient of the inverse
c        transformation, using the same schame as "coeff_f",
f        transformation, using the same schame as "COEFF_F",
*        except that "inputs" and "outputs" are transposed.
c        If "ncoeff_i" is zero, a NULL pointer may be supplied for "coeff_i".
c     lbnd_f
f     LBND_F( * ) = DOUBLE PRECISION (Given)
*        An array containing the lower bounds of the input bounding box within
*        which the ChebyMap is defined. This argument is not used or
*        accessed if
c        ncoeff_f is zero, and so a NULL pointer may be supplied.
f        NCOEFF_F is zero.
*        If supplied, the array should contain
c        "nin" elements.
f        "NIN" elements.
c     ubnd_f
f     UBND_F( * ) = DOUBLE PRECISION (Given)
*        An array containing the upper bounds of the input bounding box within
*        which the ChebyMap is defined. This argument is not used or
*        accessed if
c        ncoeff_f is zero, and so a NULL pointer may be supplied.
f        NCOEFF_F is zero.
*        If supplied, the array should contain
c        "nin" elements.
f        "NIN" elements.
c     lbnd_i
f     LBND_I( * ) = DOUBLE PRECISION (Given)
*        An array containing the lower bounds of the output bounding box within
*        which the ChebyMap is defined. This argument is not used or
*        accessed if
c        ncoeff_i is zero, and so a NULL pointer may be supplied.
f        NCOEFF_I is zero.
*        If supplied, the array should contain
c        "nout" elements.
f        "NOUT" elements.
c     ubnd_i
f     UBND_I( * ) = DOUBLE PRECISION (Given)
*        An array containing the upper bounds of the output bounding box within
*        which the ChebyMap is defined. This argument is not used or
*        accessed if
c        ncoeff_i is zero, and so a NULL pointer may be supplied.
f        NCOEFF_I is zero.
*        If supplied, the array should contain
c        "nout" elements.
f        "NOUT" elements.
c     options
f     OPTIONS = CHARACTER * ( * ) (Given)
c        Pointer to a null-terminated string containing an optional
c        comma-separated list of attribute assignments to be used for
c        initialising the new ChebyMap. The syntax used is identical to
c        that for the astSet function and may include "printf" format
c        specifiers identified by "%" symbols in the normal way.
f        A character string containing an optional comma-separated
f        list of attribute assignments to be used for initialising the
f        new ChebyMap. The syntax used is identical to that for the
f        AST_SET routine.
c     ...
c        If the "options" string contains "%" format specifiers, then
c        an optional list of additional arguments may follow it in
c        order to supply values to be substituted for these
c        specifiers. The rules for supplying these are identical to
c        those for the astSet function (and for the C "printf"
c        function).
f     STATUS = INTEGER (Given and Returned)
f        The global status.

*  Returned Value:
c     astChebyMap()
f     AST_CHEBYMAP = INTEGER
*        A pointer to the new ChebyMap.

*  Notes:
*     - A null Object pointer (AST__NULL) will be returned if this
c     function is invoked with the AST error status set, or if it
f     function is invoked with STATUS set to an error value, or if it
*     should fail for any reason.
*--
*/

/* Local Variables: */
   astDECLARE_GLOBALS          /* Pointer to thread-specific global data */
   AstChebyMap *new;            /* Pointer to new ChebyMap */
   va_list args;               /* Variable argument list */

/* Check the global status. */
   if ( !astOK ) return NULL;

/* Get a pointer to the thread specific global data structure. */
   astGET_GLOBALS(NULL);

/* Initialise the ChebyMap, allocating memory and initialising the
   virtual function table as well if necessary. */
   new = astInitChebyMap( NULL, sizeof( AstChebyMap ), !class_init,
                          &class_vtab, "ChebyMap", nin, nout,
                          ncoeff_f, coeff_f, ncoeff_i, coeff_i,
                          lbnd_f, ubnd_f, lbnd_i, ubnd_i );

/* If successful, note that the virtual function table has been
   initialised. */
   if ( astOK ) {
      class_init = 1;

/* Obtain the variable argument list and pass it along with the options string
   to the astVSet method to initialise the new ChebyMap's attributes. */
      va_start( args, status );
      astVSet( new, options, NULL, args );
      va_end( args );

/* If an error occurred, clean up by deleting the new object. */
      if ( !astOK ) new = astDelete( new );
   }

/* Return a pointer to the new ChebyMap. */
   return new;
}

AstChebyMap *astChebyMapId_( int nin, int nout, int ncoeff_f, const double coeff_f[],
                             int ncoeff_i, const double coeff_i[],
                             const double lbnd_f[], const double ubnd_f[],
                             const double lbnd_i[], const double ubnd_i[],
                             const char *options, ... ){
/*
*  Name:
*     astChebyMapId_

*  Purpose:
*     Create a ChebyMap.

*  Type:
*     Private function.

*  Synopsis:
*     #include "chebymap.h"
*     AstChebyMap *astChebyMap( int nin, int nout, int ncoeff_f, const double coeff_f[],
*                             int ncoeff_i, const double coeff_i[], const
*                             double lbnd_f[], const double ubnd_f[],
*                             double lbnd_i[], const double ubnd_i[],
*                             const char *options, ... )

*  Class Membership:
*     ChebyMap constructor.

*  Description:
*     This function implements the external (public) interface to the
*     astChebyMap constructor function. It returns an ID value (instead
*     of a true C pointer) to external users, and must be provided
*     because astChebyMap_ has a variable argument list which cannot be
*     encapsulated in a macro (where this conversion would otherwise
*     occur).
*
*     The variable argument list also prevents this function from
*     invoking astChebyMap_ directly, so it must be a re-implementation
*     of it in all respects, except for the final conversion of the
*     result to an ID value.

*  Parameters:
*     As for astChebyMap_.

*  Returned Value:
*     The ID value associated with the new ChebyMap.
*/

/* Local Variables: */
   astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
   AstChebyMap *new;              /* Pointer to new ChebyMap */
   va_list args;                 /* Variable argument list */
   int *status;                  /* Pointer to inherited status value */

/* Get a pointer to the inherited status value. */
   status = astGetStatusPtr;

/* Get a pointer to the thread specific global data structure. */
   astGET_GLOBALS(NULL);

/* Check the global status. */
   if ( !astOK ) return NULL;

/* Initialise the ChebyMap, allocating memory and initialising the
   virtual function table as well if necessary. */
   new = astInitChebyMap( NULL, sizeof( AstChebyMap ), !class_init,
                         &class_vtab, "ChebyMap", nin, nout,
                         ncoeff_f, coeff_f, ncoeff_i, coeff_i,
                         lbnd_f, ubnd_f, lbnd_i, ubnd_i );

/* If successful, note that the virtual function table has been
   initialised. */
   if ( astOK ) {
      class_init = 1;

/* Obtain the variable argument list and pass it along with the options string
   to the astVSet method to initialise the new ChebyMap's attributes. */
      va_start( args, options );
      astVSet( new, options, NULL, args );
      va_end( args );

/* If an error occurred, clean up by deleting the new object. */
      if ( !astOK ) new = astDelete( new );
   }

/* Return an ID value for the new ChebyMap. */
   return astMakeId( new );
}

AstChebyMap *astInitChebyMap_( void *mem, size_t size, int init,
                             AstChebyMapVtab *vtab, const char *name,
                             int nin, int nout, int ncoeff_f, const double coeff_f[],
                             int ncoeff_i, const double coeff_i[],
                             const double lbnd_f[], const double ubnd_f[],
                             const double lbnd_i[], const double ubnd_i[],
                             int *status ){
/*
*+
*  Name:
*     astInitChebyMap

*  Purpose:
*     Initialise a ChebyMap.

*  Type:
*     Protected function.

*  Synopsis:
*     #include "chebymap.h"
*     AstChebyMap *astInitChebyMap( void *mem, size_t size, int init,
*                                 AstChebyMapVtab *vtab, const char *name,
*                                 int nin, int nout, int ncoeff_f,
*                                 const double coeff_f[], int ncoeff_i,
*                                 const double coeff_i[]
*                                 const double lbnd_f[], const double ubnd_f[],
*                                 const double lbnd_i[], const double ubnd_i[] )

*  Class Membership:
*     ChebyMap initialiser.

*  Description:
*     This function is provided for use by class implementations to initialise
*     a new ChebyMap object. It allocates memory (if necessary) to accommodate
*     the ChebyMap plus any additional data associated with the derived class.
*     It then initialises a ChebyMap structure at the start of this memory. If
*     the "init" flag is set, it also initialises the contents of a virtual
*     function table for a ChebyMap at the start of the memory passed via the
*     "vtab" parameter.

*  Parameters:
*     mem
*        A pointer to the memory in which the ChebyMap is to be initialised.
*        This must be of sufficient size to accommodate the ChebyMap data
*        (sizeof(ChebyMap)) plus any data used by the derived class. If a value
*        of NULL is given, this function will allocate the memory itself using
*        the "size" parameter to determine its size.
*     size
*        The amount of memory used by the ChebyMap (plus derived class data).
*        This will be used to allocate memory if a value of NULL is given for
*        the "mem" parameter. This value is also stored in the ChebyMap
*        structure, so a valid value must be supplied even if not required for
*        allocating memory.
*     init
*        A logical flag indicating if the ChebyMap's virtual function table is
*        to be initialised. If this value is non-zero, the virtual function
*        table will be initialised by this function.
*     vtab
*        Pointer to the start of the virtual function table to be associated
*        with the new ChebyMap.
*     name
*        Pointer to a constant null-terminated character string which contains
*        the name of the class to which the new object belongs (it is this
*        pointer value that will subsequently be returned by the astGetClass
*        method).
*     nin
*        The number of input coordinate values per point. This is the
*        same as the number of columns in the matrix.
*     nout
*        The number of output coordinate values per point. This is the
*        same as the number of rows in the matrix.
*     ncoeff_f
*        The number of non-zero coefficients necessary to define the
*        forward transformation of the ChebyMap. If zero is supplied, the
*        forward transformation will be undefined.
*     coeff_f
*        An array containing "ncoeff_f*( 2 + nin )" elements. Each group
*	 of "2 + nin" adjacent elements describe a single coefficient of
*	 the forward transformation. Within each such group, the first
*	 element is the coefficient value; the next element is the
*	 integer index of the ChebyMap output which uses the coefficient
*	 within its defining polynomial (the first output has index 1);
*	 the remaining elements of the group give the integer powers to
*	 use with each input coordinate value (powers must not be
*	 negative)
*
*        For instance, if the ChebyMap has 3 inputs and 2 outputs, each group
*        consisting of 5 elements, A groups such as "(1.2, 2.0, 1.0, 3.0, 0.0)"
*        describes a coefficient with value 1.2 which is used within the
*        definition of output 2. The output value is incremented by the
*        product of the coefficient value, the value of input coordinate
*        1 raised to the power 1, and the value of input coordinate 2 raised
*        to the power 3. Input coordinate 3 is not used since its power is
*        specified as zero. As another example, the group "(-1.0, 1.0,
*        0.0, 0.0, 0.0 )" describes adds a constant value -1.0 onto
*        output 1 (it is a constant value since the power for every input
*        axis is given as zero).
*
*        Each final output coordinate value is the sum of the "ncoeff_f" terms
*        described by the "ncoeff_f" groups within the supplied array.
*     ncoeff_i
*        The number of non-zero coefficients necessary to define the
*        inverse transformation of the ChebyMap. If zero is supplied, the
*        inverse transformation will be undefined.
*     coeff_i
*        An array containing
*        "ncoeff_i*( 2 + nout )" elements. Each group of "2 + nout"
*        adjacent elements describe a single coefficient of the inverse
*        transformation, using the same schame as "coeff_f", except that
*        "inputs" and "outputs" are transposed.
*     lbnd_f
*        An array containing the lower bounds of the input bounding box within
*        which the ChebyMap is defined. The array should contain "nin" elements.
*     ubnd_f
*        An array containing the upper bounds of the input bounding box within
*        which the ChebyMap is defined. The array should contain "nin" elements.
*     lbnd_i
*        An array containing the lower bounds of the output bounding box within
*        which the ChebyMap is defined. The array should contain "nout" elements.
*     ubnd_i
*        An array containing the upper bounds of the output bounding box within
*        which the ChebyMap is defined. The array should contain "nout" elements.

*  Returned Value:
*     A pointer to the new ChebyMap.

*  Notes:
*     -  A null pointer will be returned if this function is invoked with the
*     global error status set, or if it should fail for any reason.
*-
*/

/* Local Variables: */
   AstChebyMap *new;
   int i;

/* Check the global status. */
   if ( !astOK ) return NULL;

/* If necessary, initialise the virtual function table. */
   if ( init ) astInitChebyMapVtab( vtab, name );

/* Initialise a PolyMap structure (the parent class) as the first component
   within the ChebyMap structure, allocating memory if necessary. */
   new = (AstChebyMap *) astInitPolyMap( mem, size, 0,
                                        (AstPolyMapVtab *) vtab, name,
                                        nin, nout, ncoeff_f, coeff_f,
                                        ncoeff_i, coeff_i );
   if ( astOK ) {

/* Initialise the ChebyMap data. */
/* ---------------------------- */

/* First initialise the pointers in case of errors. */
      new->scale_f = NULL;
      new->offset_f = NULL;
      new->scale_i = NULL;
      new->offset_i = NULL;

/* Calculate the scales and offsets that map the supplied input bounding box
   onto the range [-1,+1] on each input axis, and store them. */
      if( ncoeff_f > 0 ) {
         new->scale_f = (double *) astMalloc( sizeof( double )*nin );
         new->offset_f = (double *) astMalloc( sizeof( double )*nin );
         if( astOK ) {
            for( i = 0; i < nin; i++ ) {
               if( ubnd_f[ i ] != lbnd_f[ i ] ) {
                  new->scale_f[ i ] = 2.0/( ubnd_f[ i ] - lbnd_f[ i ] );
                  new->offset_f[ i ] = -( ubnd_f[ i ] + lbnd_f[ i ] )/( ubnd_f[ i ] - lbnd_f[ i ] );
               } else if( astOK ){
                  astError( AST__BADBX, "astInitChebyMap(%s): Input bounding box "
                            "has zero width on input axis %d.", status, name, i + 1 );
                  break;
               }
            }
         }
      }

/* Calculate the scales and offsets that map the supplied output bounding box
   onto the range [-1,+1] on each output axis, and store them. */
      if( ncoeff_i > 0 ) {
         new->scale_i = (double *) astMalloc( sizeof( double )*nout );
         new->offset_i = (double *) astMalloc( sizeof( double )*nout );
         if( astOK ) {
            for( i = 0; i < nout; i++ ) {
               if( ubnd_i[ i ] != lbnd_i[ i ] ) {
                  new->scale_i[ i ] = 2.0/( ubnd_i[ i ] - lbnd_i[ i ] );
                  new->offset_i[ i ] = -( ubnd_i[ i ] + lbnd_i[ i ] )/( ubnd_i[ i ] - lbnd_i[ i ] );
               } else if( astOK ){
                  astError( AST__BADBX, "astInitChebyMap(%s): Output bounding box "
                            "has zero width on output axis %d.", status, name, i + 1 );
                  break;
               }
            }
         }
      }

/* If an error occurred, clean up by deleting the new ChebyMap. */
      if ( !astOK ) new = astDelete( new );
   }

/* Return a pointer to the new ChebyMap. */
   return new;
}

AstChebyMap *astLoadChebyMap_( void *mem, size_t size,
                               AstChebyMapVtab *vtab, const char *name,
                               AstChannel *channel, int *status ) {
/*
*+
*  Name:
*     astLoadChebyMap

*  Purpose:
*     Load a ChebyMap.

*  Type:
*     Protected function.

*  Synopsis:
*     #include "chebymap.h"
*     AstChebyMap *astLoadChebyMap( void *mem, size_t size,
*                                   AstChebyMapVtab *vtab, const char *name,
*                                   AstChannel *channel )

*  Class Membership:
*     ChebyMap loader.

*  Description:
*     This function is provided to load a new ChebyMap using data read
*     from a Channel. It first loads the data used by the parent class
*     (which allocates memory if necessary) and then initialises a
*     ChebyMap structure in this memory, using data read from the input
*     Channel.
*
*     If the "init" flag is set, it also initialises the contents of a
*     virtual function table for a ChebyMap at the start of the memory
*     passed via the "vtab" parameter.


*  Parameters:
*     mem
*        A pointer to the memory into which the ChebyMap is to be
*        loaded.  This must be of sufficient size to accommodate the
*        ChebyMap data (sizeof(ChebyMap)) plus any data used by derived
*        classes. If a value of NULL is given, this function will
*        allocate the memory itself using the "size" parameter to
*        determine its size.
*     size
*        The amount of memory used by the ChebyMap (plus derived class
*        data).  This will be used to allocate memory if a value of
*        NULL is given for the "mem" parameter. This value is also
*        stored in the ChebyMap structure, so a valid value must be
*        supplied even if not required for allocating memory.
*
*        If the "vtab" parameter is NULL, the "size" value is ignored
*        and sizeof(AstChebyMap) is used instead.
*     vtab
*        Pointer to the start of the virtual function table to be
*        associated with the new ChebyMap. If this is NULL, a pointer
*        to the (static) virtual function table for the ChebyMap class
*        is used instead.
*     name
*        Pointer to a constant null-terminated character string which
*        contains the name of the class to which the new object
*        belongs (it is this pointer value that will subsequently be
*        returned by the astGetClass method).
*
*        If the "vtab" parameter is NULL, the "name" value is ignored
*        and a pointer to the string "ChebyMap" is used instead.

*  Returned Value:
*     A pointer to the new ChebyMap.

*  Notes:
*     - A null pointer will be returned if this function is invoked
*     with the global error status set, or if it should fail for any
*     reason.
*-
*/

#define KEY_LEN 50               /* Maximum length of a keyword */

   astDECLARE_GLOBALS            /* Pointer to thread-specific global data */
/* Local Variables: */
   AstChebyMap *new;              /* Pointer to the new ChebyMap */
   char buff[ KEY_LEN + 1 ];     /* Buffer for keyword string */
   int i;                        /* Loop index */
   int ngood;                    /* No. of non-bad values */
   int nin;                      /* No. of input coords */
   int nout;                     /* No. of output coords */

/* Get a pointer to the thread specific global data structure. */
   astGET_GLOBALS(channel);

/* Initialise. */
   new = NULL;

/* Check the global error status. */
   if ( !astOK ) return new;

/* If a NULL virtual function table has been supplied, then this is
   the first loader to be invoked for this ChebyMap. In this case the
   ChebyMap belongs to this class, so supply appropriate values to be
   passed to the parent class loader (and its parent, etc.). */
   if ( !vtab ) {
      size = sizeof( AstChebyMap );
      vtab = &class_vtab;
      name = "ChebyMap";

/* If required, initialise the virtual function table for this class. */
      if ( !class_init ) {
         astInitChebyMapVtab( vtab, name );
         class_init = 1;
      }
   }

/* Invoke the parent class loader to load data for all the ancestral
   classes of the current one, returning a pointer to the resulting
   partly-built ChebyMap. */
   new = astLoadPolyMap( mem, size, (AstPolyMapVtab *) vtab, name,
                         channel );

   if ( astOK ) {

/* Get the number of inputs and outputs for the uninverted Mapping. */
   nin = ( (AstMapping *) new )->nin;
   nout = ( (AstMapping *) new )->nout;

/* Initialise values */
   new->scale_f = NULL;
   new->offset_f = NULL;
   new->scale_i = NULL;
   new->offset_i = NULL;

/* Read input data. */
/* ================ */

/* Request the input Channel to read all the input data appropriate to
   this class into the internal "values list". */
      astReadClassData( channel, "ChebyMap" );

/* Is the forward transformation defined? */
      if( ((AstPolyMap *) new)->ncoeff_f ) {

/* Allocate memory to hold the scales and offsets for the input axes. */
         new->scale_f = astMalloc( sizeof( double )*(size_t) nin );
         new->offset_f = astMalloc( sizeof( double )*(size_t) nin );
         if( astOK ) {

/* Get the scale factors. */
            ngood = 0;
            for( i = 0; i < nin; i++ ){
               (void) sprintf( buff, "fscl%d", i + 1 );
               (new->scale_f)[ i ] = astReadDouble( channel, buff, AST__BAD );
               if( (new->scale_f)[ i ] != AST__BAD ) ngood++;
            }

/* Get the offsets of the bounding box. */
            for( i = 0; i < nin; i++ ){
               (void) sprintf( buff, "foff%d", i + 1 );
               (new->offset_f)[ i ] = astReadDouble( channel, buff, AST__BAD );
               if( (new->offset_f)[ i ] != AST__BAD ) ngood++;
            }

/* The scale and offset values should all be AST__BAD if the transformation
   is a standard polynomial. Anull the scale and offset arrays to
   indicate this. */
            if( ngood == 0 ) {
               new->scale_f = astFree( new->scale_f );
               new->offset_f = astFree( new->offset_f );

/* Otherwise, there should be no bad values. */
            } else if( ngood != 2*nin && astOK ) {
               astError( AST__OBJIN, "astLoadChebyMap: insufficient scale "
                         "and offset values for the forward transformation "
                         "in loaded ChebyMap.", status );
            }
         }
      }

/* Is the inverse transformation defined? */
      if( ((AstPolyMap *) new)->ncoeff_i ) {

/* Allocate memory to hold the scales and offsets for the output axes. */
         new->scale_i = astMalloc( sizeof( double )*(size_t) nout );
         new->offset_i = astMalloc( sizeof( double )*(size_t) nout );
         if( astOK ) {

/* Get the scale factors. */
            ngood = 0;
            for( i = 0; i < nout; i++ ){
               (void) sprintf( buff, "iscl%d", i + 1 );
               (new->scale_i)[ i ] = astReadDouble( channel, buff, AST__BAD );
               if( (new->scale_i)[ i ] != AST__BAD ) ngood++;
            }

/* Get the offsets of the bounding box. */
            for( i = 0; i < nout; i++ ){
               (void) sprintf( buff, "ioff%d", i + 1 );
               (new->offset_i)[ i ] = astReadDouble( channel, buff, AST__BAD );
               if( (new->offset_i)[ i ] != AST__BAD ) ngood++;
            }

/* The scale and offset values should all be AST__BAD if the transformation
   is a standard polynomial. Anull the scale and offset arrays to
   indicate this. */
            if( ngood == 0 ) {
               new->scale_i = astFree( new->scale_i );
               new->offset_i = astFree( new->offset_i );

/* Otherwise, there should be no bad values. */
            } else if( ngood != 2*nout && astOK ) {
               astError( AST__OBJIN, "astLoadChebyMap: insufficient scale "
                         "and offset values for the inverse transformation "
                         "in loaded ChebyMap.", status );
            }
         }
      }

/* If an error occurred, clean up by deleting the new ChebyMap. */
      if ( !astOK ) new = astDelete( new );
   }

/* Return the new ChebyMap pointer. */
   return new;

/* Undefine macros local to this function. */
#undef KEY_LEN
}

/* Virtual function interfaces. */
/* ============================ */
/* These provide the external interface to the virtual functions defined by
   this class. Each simply checks the global error status and then locates and
   executes the appropriate member function, using the function pointer stored
   in the object's virtual function table (this pointer is located using the
   astMEMBER macro defined in "object.h").

   Note that the member function may not be the one defined here, as it may
   have been over-ridden by a derived class. However, it should still have the
   same interface. */


void astChebyDomain_( AstChebyMap *this, int forward, double *lbnd, double *ubnd, int *status ){
   if ( !astOK ) return;
   (**astMEMBER(this,ChebyMap,ChebyDomain))( this, forward, lbnd, ubnd, status );
}