File: uves_utils_wrappers.c

package info (click to toggle)
cpl-plugin-uves 6.1.3+dfsg-2
  • links: PTS, VCS
  • area: main
  • in suites: bullseye, sid
  • size: 23,128 kB
  • sloc: ansic: 171,056; sh: 4,359; python: 3,002; makefile: 1,322
file content (3076 lines) | stat: -rw-r--r-- 104,775 bytes parent folder | download
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
2405
2406
2407
2408
2409
2410
2411
2412
2413
2414
2415
2416
2417
2418
2419
2420
2421
2422
2423
2424
2425
2426
2427
2428
2429
2430
2431
2432
2433
2434
2435
2436
2437
2438
2439
2440
2441
2442
2443
2444
2445
2446
2447
2448
2449
2450
2451
2452
2453
2454
2455
2456
2457
2458
2459
2460
2461
2462
2463
2464
2465
2466
2467
2468
2469
2470
2471
2472
2473
2474
2475
2476
2477
2478
2479
2480
2481
2482
2483
2484
2485
2486
2487
2488
2489
2490
2491
2492
2493
2494
2495
2496
2497
2498
2499
2500
2501
2502
2503
2504
2505
2506
2507
2508
2509
2510
2511
2512
2513
2514
2515
2516
2517
2518
2519
2520
2521
2522
2523
2524
2525
2526
2527
2528
2529
2530
2531
2532
2533
2534
2535
2536
2537
2538
2539
2540
2541
2542
2543
2544
2545
2546
2547
2548
2549
2550
2551
2552
2553
2554
2555
2556
2557
2558
2559
2560
2561
2562
2563
2564
2565
2566
2567
2568
2569
2570
2571
2572
2573
2574
2575
2576
2577
2578
2579
2580
2581
2582
2583
2584
2585
2586
2587
2588
2589
2590
2591
2592
2593
2594
2595
2596
2597
2598
2599
2600
2601
2602
2603
2604
2605
2606
2607
2608
2609
2610
2611
2612
2613
2614
2615
2616
2617
2618
2619
2620
2621
2622
2623
2624
2625
2626
2627
2628
2629
2630
2631
2632
2633
2634
2635
2636
2637
2638
2639
2640
2641
2642
2643
2644
2645
2646
2647
2648
2649
2650
2651
2652
2653
2654
2655
2656
2657
2658
2659
2660
2661
2662
2663
2664
2665
2666
2667
2668
2669
2670
2671
2672
2673
2674
2675
2676
2677
2678
2679
2680
2681
2682
2683
2684
2685
2686
2687
2688
2689
2690
2691
2692
2693
2694
2695
2696
2697
2698
2699
2700
2701
2702
2703
2704
2705
2706
2707
2708
2709
2710
2711
2712
2713
2714
2715
2716
2717
2718
2719
2720
2721
2722
2723
2724
2725
2726
2727
2728
2729
2730
2731
2732
2733
2734
2735
2736
2737
2738
2739
2740
2741
2742
2743
2744
2745
2746
2747
2748
2749
2750
2751
2752
2753
2754
2755
2756
2757
2758
2759
2760
2761
2762
2763
2764
2765
2766
2767
2768
2769
2770
2771
2772
2773
2774
2775
2776
2777
2778
2779
2780
2781
2782
2783
2784
2785
2786
2787
2788
2789
2790
2791
2792
2793
2794
2795
2796
2797
2798
2799
2800
2801
2802
2803
2804
2805
2806
2807
2808
2809
2810
2811
2812
2813
2814
2815
2816
2817
2818
2819
2820
2821
2822
2823
2824
2825
2826
2827
2828
2829
2830
2831
2832
2833
2834
2835
2836
2837
2838
2839
2840
2841
2842
2843
2844
2845
2846
2847
2848
2849
2850
2851
2852
2853
2854
2855
2856
2857
2858
2859
2860
2861
2862
2863
2864
2865
2866
2867
2868
2869
2870
2871
2872
2873
2874
2875
2876
2877
2878
2879
2880
2881
2882
2883
2884
2885
2886
2887
2888
2889
2890
2891
2892
2893
2894
2895
2896
2897
2898
2899
2900
2901
2902
2903
2904
2905
2906
2907
2908
2909
2910
2911
2912
2913
2914
2915
2916
2917
2918
2919
2920
2921
2922
2923
2924
2925
2926
2927
2928
2929
2930
2931
2932
2933
2934
2935
2936
2937
2938
2939
2940
2941
2942
2943
2944
2945
2946
2947
2948
2949
2950
2951
2952
2953
2954
2955
2956
2957
2958
2959
2960
2961
2962
2963
2964
2965
2966
2967
2968
2969
2970
2971
2972
2973
2974
2975
2976
2977
2978
2979
2980
2981
2982
2983
2984
2985
2986
2987
2988
2989
2990
2991
2992
2993
2994
2995
2996
2997
2998
2999
3000
3001
3002
3003
3004
3005
3006
3007
3008
3009
3010
3011
3012
3013
3014
3015
3016
3017
3018
3019
3020
3021
3022
3023
3024
3025
3026
3027
3028
3029
3030
3031
3032
3033
3034
3035
3036
3037
3038
3039
3040
3041
3042
3043
3044
3045
3046
3047
3048
3049
3050
3051
3052
3053
3054
3055
3056
3057
3058
3059
3060
3061
3062
3063
3064
3065
3066
3067
3068
3069
3070
3071
3072
3073
3074
3075
3076
/*                                                                              *
 *   This file is part of the ESO UVES Pipeline                                 *
 *   Copyright (C) 2004,2005 European Southern Observatory                      *
 *                                                                              *
 *   This library is free software; you can redistribute it and/or modify       *
 *   it under the terms of the GNU General Public License as published by       *
 *   the Free Software Foundation; either version 2 of the License, or          *
 *   (at your option) any later version.                                        *
 *                                                                              *
 *   This program is distributed in the hope that it will be useful,            *
 *   but WITHOUT ANY WARRANTY; without even the implied warranty of             *
 *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the              *
 *   GNU General Public License for more details.                               *
 *                                                                              *
 *   You should have received a copy of the GNU General Public License          *
 *   along with this program; if not, write to the Free Software                *
 *   Foundation, 51 Franklin St, Fifth Floor, Boston, MA  02111-1307  USA       *
 *                                                                              */

/*
 * $Author: amodigli $
 * $Date: 2012-05-01 06:27:56 $
 * $Revision: 1.123 $
 * $Name: not supported by cvs2svn $
 */

#ifdef HAVE_CONFIG_H
#  include <config.h>
#endif

/*----------------------------------------------------------------------------*/
/**
   @defgroup uves_utils_wrappers  Utility functions (wrappers)

   This module contains wrapper functions, convenience functions and
   simple extensions of CPL functions.
 
 */
/*----------------------------------------------------------------------------*/

#include <uves_utils_wrappers.h>

#include <uves_utils.h>
#include <uves_utils_cpl.h>
#include <uves_error.h>
#include <uves_dump.h>
#include <cpl.h>

#include <irplib_utils.h>
#include <stdbool.h>
#include <assert.h>

/*-----------------------------------------------------------------------------
                                   Local functions
 -----------------------------------------------------------------------------*/

static int get_candidate(const double *a, const int ia[],
             int M, int N, int D,
             double lambda,
             int    (*f)(const double x[], const double a[], 
                     double *result),
             int (*dfda)(const double x[], const double a[], 
                     double result[]),
             const double *x,
             const double *y,
             const double *sigma,
             double *partials,
             cpl_matrix *alpha,
             cpl_matrix *beta,
             double *a_da);

static double get_chisq(int N, int D,
            int (*f)(const double x[], const double a[], 
                 double *result),
            const double *a,
            const double *x,
            const double *y,
            const double *sigma);


static cpl_image * 
uves_image_filter_wrapper(const cpl_image *b, 
                          const cpl_matrix *k, 
                          cpl_filter_mode mode);
cpl_image * 
uves_image_filter_median(const cpl_image * img, const cpl_matrix * mx)
{
	return uves_image_filter_wrapper(img, mx, CPL_FILTER_MEDIAN);
}

cpl_image * 
uves_image_filter_linear(const cpl_image *img, const cpl_matrix * mx)
{
	return uves_image_filter_wrapper(img, mx, CPL_FILTER_LINEAR);

}

/*----------------------------------------------------------------------------
                                   Defines
 ---------------------------------------------------------------------------*/

/** Provide nicer syntax */
#define image_is_rejected(badmap, x, y) \
  ((badmap) != NULL && (badmap)[((x)-1) + ((y)-1)*nx] == CPL_BINARY_1)

#ifndef UVES_FIT_MAXITER
#define UVES_FIT_MAXITER 1000
#endif

/*-----------------------------------------------------------------------------
                                   Types
 -----------------------------------------------------------------------------*/
/* @cond */
typedef struct {
    double x;
    double y;
} uves_fit_1d_input;
/* @endcond */

/*-----------------------------------------------------------------------------
                                   Implementation
 -----------------------------------------------------------------------------*/
/**@{*/



static cpl_image * 
uves_image_filter_wrapper(const cpl_image *b, 
                          const cpl_matrix *k, 
                          cpl_filter_mode mode)
{
	const double EPSILON = 1E-5;
	int nx   = cpl_image_get_size_x(b);
	int ny   = cpl_image_get_size_y(b);
	int nrow = cpl_matrix_get_nrow(k);
	int ncol = cpl_matrix_get_ncol(k);
	int i, j;
	cpl_type type = cpl_image_get_type(b);
	cpl_image * a = cpl_image_new(nx, ny, type);
	// where m is a cpl_mask with a CPL_BINARY_1 whereever k has a 1.0.
	cpl_mask* m = cpl_mask_new(ncol, nrow);
	//cpl_msg_warning(cpl_func, "nx[%d], ny[%d], ncol[%d], nrow[%d]", nx, ny, ncol, nrow);
	for (i = 0; i < ncol ; i++)
	{
		for (j = 0; j < nrow ; j++)
		{
			double value = cpl_matrix_get(k, j, i);
			if (fabs(value - 1.0) < EPSILON)
			{
				cpl_mask_set(m, i + 1, j + 1, CPL_BINARY_1);
			}
		}
	}

	cpl_image_filter_mask(a, b, m, mode, CPL_BORDER_FILTER);
	cpl_mask_delete(m);
	return a;
}

cpl_image*
uves_image_filter_mode(const cpl_image* b,
                      const cpl_matrix * ker,
                      cpl_filter_mode filter)
{
  int nx   = cpl_image_get_size_x(b);
  int ny   = cpl_image_get_size_y(b);
  int type = cpl_image_get_type(b);
  cpl_image * a = cpl_image_new(nx, ny, type);

  switch(filter) {
  case CPL_FILTER_MEDIAN:
    check_nomsg(cpl_image_filter(a, b, ker, CPL_FILTER_MEDIAN, CPL_BORDER_FILTER));
    break;
  case CPL_FILTER_LINEAR:
    check_nomsg(cpl_image_filter(a, b, ker, CPL_FILTER_LINEAR, CPL_BORDER_FILTER));
    break;
  case CPL_FILTER_STDEV:
    cpl_image_filter(a, b, ker, CPL_FILTER_STDEV, CPL_BORDER_FILTER);
    break;
  case CPL_FILTER_MORPHO:
    cpl_image_filter(a, b, ker, CPL_FILTER_MORPHO, CPL_BORDER_FILTER);
    break;
  default:
    uves_msg_error("Filter type not supported");
    return NULL;
  }
 cleanup:

  return a;

}


/*----------------------------------------------------------------------------*/
/**
   @brief   Reject all pixels in an image
   @param   image      operand

   Note: calls accessor function on every pixels (i.e. could be optimized)
*/
/*----------------------------------------------------------------------------*/
void uves_image_reject_all(cpl_image *image)
{
    int nx, ny;
    int x, y;

    assure_nomsg( image != NULL, CPL_ERROR_NULL_INPUT );
    nx = cpl_image_get_size_x(image);
    ny = cpl_image_get_size_y(image);

    for (y = 1; y <= ny; y++) {
        for (x = 1; x <= nx; x++) {
            cpl_image_reject(image, x, y);
        }
    }
    
  cleanup:
    return;    
}


/*----------------------------------------------------------------------------*/
/**
   @brief   Apply a 1d fit to an image subrow or subcolumn
   @param   image       The image to fit.
   @param   noise       The uncertainty (one sigma, gaussian errors assumed)
                        associated with the image.
                        If NULL, constant uncertainties are assumed.
   @param   image_badmap Pointer to image bad pixel map, or NULL (all good pixels)
                         . Necessary for efficiency with CPL2
   @param   horizontal  If true, a subrow is fitted. If false, a subcolumn is
                        fitted.
   @param   fix_back    If true, the background is not fitted but kept constant
   @param   fit_back    If true, only the background and area are fitted. In
                        this case the normalization can be positive or negative,
                        and therefore a different centering method (which does
                        not assume that the normalization is positive) is used.
   @param   y_0         Row/column to fit (FITS convention).
   @param   xlo         First column/row (inclusive) to fit (FITS convention).
   @param   xhi         Last column/row (inclusive) to fit (FITS convention).
   @param   x0          (output) Center of best fit gaussian.
   @param   sigma       (output) Width of best fit gaussian. A positive number
                        on success.
   @param   norm        (output) Area of gaussian. Positive on success
   @param   background  (output) Fitted background level.
   @param   slope       (output) Fitted slope (only if M = 5)
   @param   mse         (output) If non-NULL, the mean squared error of the best
                        fit is returned.
   @param   red_chisq   (output) If non-NULL, the reduced chi square of the best
                        fit is returned. This requires the noise image to be
            specified.
   @param   covariance  (output) If non-NULL, the formal covariance matrix of
                        the best fit is returned. This requires the noise image
            to be specified. The order of fit parameters in the
            covariance matrix is defined as (@em x0, @em sigma,
            @em norm, @em background), for example the (3,3)
            element of the matrix (counting from zero) is the
            variance of the fitted @em background. The matrix must
            be deallocated by calling @c cpl_matrix_delete() . On
            error, NULL is returned.
   @param   f           The fit function
   @param   dfda        The derivative of the fit function
   @param   M           Number of fit parameters, must be 4 or 5

   @return  CPL_ERROR_NONE iff okay

   This function fits part of an image row or column to a 1d gaussian function

   If the parameter @em horizontal is set to true, the gaussian is fitted to the
   image subrow (@em xlo, @em y_0) - (@em xhi, @em y_0).
   Otherwise the subcolumn (@em y_0, @em xlo) - (@em y_0, @em xhi) is used.

   The values to fit are read from the input @em image. Optionally, a @em noise
   image (of same size as the input image) may be specified. Both images must
   have type CPL_TYPE_DOUBLE. If a pixel is bad in either the input image or in
   the noise image, it is excluded from the fit. To perform the fit, at least
   one good pixel must exist in the specified subwindow (but it is of course
   recommendable to use as wide a fitting window as possible to get a
   well-defined fit).
*/
/*----------------------------------------------------------------------------*/

cpl_error_code
uves_fit_1d_image(const cpl_image *image, const cpl_image *noise,
          const cpl_binary *image_badmap,
          bool horizontal, bool fix_back, bool fit_back,
          int xlo, int xhi, int y_0,
          double *x0, double *sigma, double *norm, double *background,
                  double *slope,
          double *mse, double *red_chisq,
          cpl_matrix **covariance,
          int (*f)   (const double x[], const double a[], double *result),
          int (*dfda)(const double x[], const double a[], double result[]),
          int M)
{
    cpl_vector *x = NULL;
    cpl_vector *y = NULL;
    cpl_vector *sigma_y = NULL;     /* Noise vector          */
    
    cpl_fit_mode fit_pattern;
    int nx, ny;                     /* Image size            */
    int N;                          /* Number of good pixels */
    int i, j;
    cpl_type image_type;

    const double *image_data       = NULL;        /* Pointer to data */
    const double *noise_data       = NULL;        /* Pointer to data */

    assure( x0 != NULL        , CPL_ERROR_NULL_INPUT, "Null fit parameter");
    assure( sigma != NULL     , CPL_ERROR_NULL_INPUT, "Null fit parameter");
    assure( norm != NULL      , CPL_ERROR_NULL_INPUT, "Null fit parameter");
    assure( background != NULL, CPL_ERROR_NULL_INPUT, "Null fit parameter");
    /* mse, red_chisq, covariance may be NULL */

    assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");
   
    image_type = cpl_image_get_type(image);

    /* To support the following types, use cpl_image_get() or
       more multiple type pointers to data buffer.
       cpl_ensure_code( image_type == CPL_TYPE_INT ||
       image_type == CPL_TYPE_FLOAT ||
       image_type == CPL_TYPE_DOUBLE,
       CPL_ERROR_TYPE_MISMATCH);
    */

    assure( image_type == CPL_TYPE_DOUBLE, CPL_ERROR_UNSUPPORTED_MODE,
        "Unsupported type: %s", uves_tostring_cpl_type(image_type));

    image_data = cpl_image_get_data_double_const(image);

    if (noise != NULL)
    {
        image_type = cpl_image_get_type(noise);
        /*  See comment above.
          cpl_ensure_code( image_type == CPL_TYPE_INT ||
          image_type == CPL_TYPE_FLOAT ||
          image_type == CPL_TYPE_DOUBLE,
          CPL_ERROR_TYPE_MISMATCH);
        */
        assure( image_type == CPL_TYPE_DOUBLE, CPL_ERROR_UNSUPPORTED_MODE, 
            "Unsupported type: %s", uves_tostring_cpl_type(image_type));

        noise_data = cpl_image_get_data_double_const(noise);
    }   

    nx = cpl_image_get_size_x(image);
    ny = cpl_image_get_size_y(image);

    if (horizontal)
    {
        assure( 1 <= xlo && xlo <= xhi && xhi <= nx &&
            1 <= y_0  && y_0  <= ny, CPL_ERROR_ACCESS_OUT_OF_RANGE,
            "Illegal window (%d, %d)-(%d, %d), image: %dx%d",
            xlo, y_0, xhi, y_0,
            nx, ny);
    }
    else
    {
        assure( 1 <= xlo && xlo <= xhi && xhi <= ny &&
            1 <= y_0  && y_0  <= nx,
            CPL_ERROR_ACCESS_OUT_OF_RANGE,
            "Illegal window (%d, %d)-(%d, %d), image: %dx%d",
            y_0, xlo, y_0, xhi,
            nx, ny);
    }
    
    /* Noise image must have same size
     * as the input image
     */
    if (noise != NULL)
    {
        assure( cpl_image_get_size_x(noise) == nx &&
            cpl_image_get_size_y(noise) == ny,
            CPL_ERROR_INCOMPATIBLE_INPUT, "Noise image: %" CPL_SIZE_FORMAT "x%" CPL_SIZE_FORMAT ", image: %dx%d:",
            cpl_image_get_size_x(noise),
            cpl_image_get_size_y(noise),
            nx, ny);
    }
   
    /* Count good pixels in sub-window, check that noise image is positive */
    N = 0;
    for (i = xlo; i <= xhi; i++)
    {
        if ( !image_is_rejected(image_badmap,
                    (horizontal) ? i : y_0,
                    (horizontal) ? y_0 : i) )
        {
            if ( noise != NULL)
            {
                if( !image_is_rejected(image_badmap,
                           (horizontal) ? i : y_0,
                           (horizontal) ? y_0 : i) )
                {
                    /* Noise image must be positive (only check
                       pixels that are actually used) */
                    assure( /*cpl_image_get(noise,
                          (horizontal) ? i : y_0,
                          (horizontal) ? y_0 : i,
                          &pis_rejected)*/
                    noise_data[(((horizontal) ? i : y_0) - 1) +
                           (((horizontal) ? y_0 : i) - 1) * nx]
                    > 0,
                    CPL_ERROR_ILLEGAL_INPUT,
                    "Non-positive noise at (%d, %d): %e",
                    (horizontal) ? i : y_0,
                    (horizontal) ? y_0 : i,
                    noise_data[(((horizontal) ? i : y_0) - 1) +
                           (((horizontal) ? y_0 : i) - 1) * nx]
                    );
                    
                    N += 1;
                }
                else
                {
                    /* Pixel value is good, but noise pixel is
                       bad. Don't use. */
                }
            }
            else
            {
                /* Pixel is good. No noise image */
                N += 1;
            }
        }
    }
   
    /* Check that there is at least one good pixel */
    assure( N >= 1, CPL_ERROR_ILLEGAL_INPUT, "Only %d good pixel(s)", N);

    /* Allocate space */
    x = cpl_vector_new(N);
    y = cpl_vector_new(N);
    if (noise != NULL)
    {
        sigma_y = cpl_vector_new(N);
        assure_mem( sigma_y );
    }

    if (fix_back)
    {
        fit_pattern = CPL_FIT_CENTROID | CPL_FIT_STDEV | CPL_FIT_AREA;
    }
    else if (fit_back)
    {
        fit_pattern = CPL_FIT_AREA | CPL_FIT_OFFSET;
    }
    else
    {
        fit_pattern = CPL_FIT_ALL;
    }
   
    assure_mem( x );
    assure_mem( y );
    
    /* Copy the N good pixels from the input image to vectors,
       j count good pixels */
    for (i = xlo, j = 0;
     i <= xhi;
     i++)
    {
        double flux;
        
        /*
          flux = cpl_image_get(image,
          (horizontal) ? xlo+i : y_0,
          (horizontal) ? y_0 : xlo+i,
          &pis_rejected);
        */
        
        flux = image_data[(((horizontal) ? i : y_0) - 1) +
                  (((horizontal) ? y_0 : i) - 1) * nx];
       
        //if (!pis_rejected)
        if ( !image_is_rejected(image_badmap,
                    (horizontal) ? i : y_0,
                    (horizontal) ? y_0 : i) )
        {
            if (noise != NULL)
            {
                double flux_noise;
               
                /* flux_noise = cpl_image_get(noise,
                   (horizontal) ? xlo+i : y_0,
                   (horizontal) ? y_0 : xlo+i,
                   &pis_rejected);
                */
               
                flux_noise = noise_data
                [(((horizontal) ? i : y_0) - 1) +
                 (((horizontal) ? y_0 : i) - 1) * nx];
               
                //if (!pis_rejected)
                if ( !image_is_rejected(image_badmap,
                            (horizontal) ?
                            i : y_0,
                            (horizontal)
                            ? y_0 : i) )
                {
                    cpl_vector_set(x,       j, i);
                    cpl_vector_set(y,       j, flux);
                    cpl_vector_set(sigma_y, j, flux_noise);
                    j++;
                }
            }
            else
            {
                cpl_vector_set(x, j, i);
                cpl_vector_set(y, j, flux);
                j++;
            }
        }
    }
    passure( j == N, "%d %d", j, N);
    
    check( uves_fit_1d(x, NULL,      /* x, sigma_x */
               y, sigma_y,
               fit_pattern, fit_back,
               x0, sigma, norm, background,
                       slope,
               mse, red_chisq,
               covariance,
               f, dfda, M),
       "Fit failed");

    
  cleanup:
    uves_free_vector(&x);
    uves_free_vector(&y);
    uves_free_vector(&sigma_y);

    return cpl_error_get_code();
}

/*----------------------------------------------------------------------------*/
/**
   @brief   Define order of input data
   @param   left      Left operand
   @param   right     Right operand

   This function uses void pointers because it is used as input for ANSI's qsort().
*/
/*----------------------------------------------------------------------------*/
static int uves_fit_1d_compare(const void *left,
                   const void *right)
{
    return 
    (((uves_fit_1d_input *)left )->x <  
     ((uves_fit_1d_input *)right)->x) ? -1 :
    (((uves_fit_1d_input *)left )->x == 
     ((uves_fit_1d_input *)right)->x) ? 0  : 1;
}

/*----------------------------------------------------------------------------*/
/** 
   @cond

   This function is just like the CPL gaussian fit function,
   except it parametrizes the function to fit.

   The number of parameters, M, must be 4 or 5.

   The meaning of fit parameters must be:
   a[0]: centroid
   a[1]: stdev
   a[2]: area
   a[3]: offset (sky)
   (a[4]: linear sky term)

   Also refer to the documentation for the CPL function.

   If @em fit_back is true, the normalization (area) is allowed
   to be positive or negative, and therefore the initial centering
   method is different. In this case @em fit_pars must
   be == FIT_OFFSET | FIT_AREA.
*/
/*----------------------------------------------------------------------------*/
cpl_error_code
uves_fit_1d(cpl_vector *x, const cpl_vector *sigma_x,
            cpl_vector *y, const cpl_vector *sigma_y,
            cpl_fit_mode fit_pars, bool fit_back,
            double *x0, double *sigma, double *area, double *offset, double *slope,
            double *mse, double *red_chisq,
            cpl_matrix **covariance,
            int (*f)   (const double x[], const double a[], double *result),
            int (*dfda)(const double x[], const double a[], double result[]),
            int M)
{
    cpl_matrix *x_matrix = NULL; /* LM algorithm needs a matrix,
                      not a vector                 */

    int N;                          /* Number of data points        */


    /* Initial parameter values */
    double x0_guess    = 0;  /* Avoid warnings about uninitialized variables */
    double sigma_guess = 0;
    double area_guess;
    double offset_guess;

    cpl_ensure_code( M == 4 || M == 5, CPL_ERROR_UNSUPPORTED_MODE);

    /* Validate input */
    cpl_ensure_code( x       != NULL, CPL_ERROR_NULL_INPUT);
    cpl_ensure_code( sigma_x == NULL, CPL_ERROR_UNSUPPORTED_MODE);
    cpl_ensure_code( y       != NULL, CPL_ERROR_NULL_INPUT);
    /* sigma_y may be NULL or non-NULL */
    
    N = cpl_vector_get_size(x);

    cpl_ensure_code( N == cpl_vector_get_size(y),
             CPL_ERROR_INCOMPATIBLE_INPUT);

    if (sigma_x != NULL)
    {
        cpl_ensure_code( N == cpl_vector_get_size(sigma_x),
                 CPL_ERROR_INCOMPATIBLE_INPUT);
    }
    if (sigma_y != NULL)
    {
        cpl_ensure_code( N == cpl_vector_get_size(sigma_y),
                 CPL_ERROR_INCOMPATIBLE_INPUT);
    }

    cpl_ensure_code( x0     != NULL &&
             sigma  != NULL &&
             area   != NULL &&
             offset != NULL &&
                     (M != 5 || slope != NULL), CPL_ERROR_NULL_INPUT);

    if (! (fit_pars & CPL_FIT_STDEV))
    {
        cpl_ensure_code( *sigma > 0, CPL_ERROR_ILLEGAL_INPUT);
    }

    cpl_ensure_code( !fit_back || fit_pars == (CPL_FIT_OFFSET | CPL_FIT_AREA),
             CPL_ERROR_INCOMPATIBLE_INPUT);

    /* Input area must be positive if fit_back is false */
    if (! (fit_pars & CPL_FIT_AREA) && !fit_back)
    {
        cpl_ensure_code( *area > 0, CPL_ERROR_ILLEGAL_INPUT);
    }

    /* mse, red_chisq may be NULL */

    /* Need more than number_of_parameters points to calculate chi^2.
     * There are less than 5 parameters. */
    cpl_ensure_code( red_chisq == NULL || N >= 5, CPL_ERROR_ILLEGAL_INPUT);
    
    if (covariance != NULL) *covariance = NULL;
    /* If covariance computation is requested, then
     * return either the covariance matrix or NULL
     * (don't return undefined pointer).
     */
    
    /* Cannot compute chi square & covariance without sigma_y */
    cpl_ensure_code( (red_chisq == NULL && covariance == NULL) || 
             sigma_y != NULL,
             CPL_ERROR_INCOMPATIBLE_INPUT);
    
    /* Create matrix from x-data */
    x_matrix = cpl_matrix_wrap(N, 1, cpl_vector_get_data(x));
    if (x_matrix == NULL)
    {
        cpl_ensure_code(
                CPL_FALSE,
        CPL_ERROR_ILLEGAL_OUTPUT);
    }
    
    /* Check that any provided sigmas are positive. */
    if (sigma_x != NULL && cpl_vector_get_min(sigma_x) <= 0)
    {
        cpl_matrix_unwrap(x_matrix);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_INPUT);
    }
    if (sigma_y != NULL && cpl_vector_get_min(sigma_y) <= 0)
    {
        cpl_matrix_unwrap(x_matrix);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_INPUT);
    }

    /* Compute guess parameters using robust estimation
     * (This step might be improved by taking into account the 
     * uncertainties but for simplicity's sake do not)
     */
    if (fit_back)
    {
        /* We need to estimate only these two parameters */
        assert( fit_pars == CPL_FIT_OFFSET || CPL_FIT_AREA);

        offset_guess = cpl_vector_get_median_const(y);
        area_guess = N*(cpl_vector_get_mean(y) - offset_guess);
        /* Sum of (flux-offset) */

        x0_guess = *x0;
        sigma_guess = *sigma;
    }
    else {
    double sum  = 0.0;
    double quartile[3];
    double fraction[3] = {0.25 , 0.50 , 0.75};
    const double *y_data = cpl_vector_get_data_const(y);

    if (fit_pars & CPL_FIT_OFFSET)
        {
        /* Estimate offset as 25% percentile of y-values.
           (The minimum value may be too low for noisy input,
           the median is too high if there is not much
           background in the supplied data, so use
           something inbetween).
        */

        cpl_vector *y_dup = cpl_vector_duplicate(y);
        
        if (y_dup == NULL)
            {
            cpl_matrix_unwrap(x_matrix);
            cpl_ensure_code(
                CPL_FALSE,
                CPL_ERROR_ILLEGAL_OUTPUT);
            }
        
        offset_guess = uves_utils_get_kth_double(
            cpl_vector_get_data(y_dup), N, N/4);
        
        cpl_vector_delete(y_dup);
        }
    else
        {
        offset_guess = *offset;
        }
    
    /* Get quartiles of distribution
       (only bother if it's needed for estimation of x0 or sigma) */
    if ( (fit_pars & CPL_FIT_CENTROID) ||
         (fit_pars & CPL_FIT_STDEV   )
        )
        {
        /* The algorithm requires the input to be sorted
           as function of x, so do that (using qsort), and
           work on the sorted copy. Of course, the y-vector
           must be re-ordered along with x.
           sigma_x and sigma_y are not used, so don't copy those.
        */
        
        uves_fit_1d_input
            *sorted_input = cpl_malloc(N * sizeof(*sorted_input));
        const double *x_data = cpl_matrix_get_data_const(x_matrix);
        cpl_boolean is_sorted = CPL_TRUE;
        int i;
        
        if (sorted_input == NULL)
            {
            cpl_matrix_unwrap(x_matrix);
            cpl_ensure_code(
                CPL_FALSE,
                CPL_ERROR_ILLEGAL_INPUT);
            }
        
        for (i = 0; i < N; i++)
            {
            sorted_input[i].x = x_data[i];
            sorted_input[i].y = y_data[i];
            
            is_sorted = is_sorted && 
                (i==0 || (x_data[i-1] < x_data[i]));
            }            
        
        if (!is_sorted)
            {
            qsort(sorted_input, N, sizeof(*sorted_input),
                  &uves_fit_1d_compare);
            }

        for(i = 0; i < N; i++)
            {
            double flux = sorted_input[i].y;
            
            sum += (flux - offset_guess);
            }
        /* Note that 'sum' must be calculated the same way as
           'running_sum' below, Otherwise (due to round-off error)
           'running_sum' might end up being different from 'sum'(!).
           Specifically, it will not work to calculate 'sum' as
           
           (flux1 + ... + fluxN)  -  N*offset_guess
        */
        
        if (sum > 0.0)
            {
            double flux, x1, x2;
            double running_sum = 0.0;
            int j;
            
            i = 0;
            flux = sorted_input[i].y - offset_guess;
            
            for (j = 0; j < 3; j++)
                {
                double limit = fraction[j] * sum;
                double k;
                
                while (running_sum + flux < limit && i < N-1)
                    {
                    running_sum += flux;
                    i++;
                    flux = sorted_input[i].y - offset_guess;
                    }

                /* Fraction [0;1] of current flux needed
                   to reach the quartile */
                k = (limit - running_sum)/flux;
                
                if (k <= 0.5)
                    {
                    /* Interpolate linearly between
                       current and previous position
                    */
                    if (0 < i)
                        {
                        x1 = sorted_input[i-1].x;
                        x2 = sorted_input[i].x;
                        
                        quartile[j] = 
                            x1*(0.5-k) +
                            x2*(0.5+k);
                        /*
                          k=0   => quartile = midpoint,
                          k=0.5 => quartile = x2
                        */
                        }
                    else
                        {
                        quartile[j] = sorted_input[i].x;
                        }
                    }
                else
                    {
                    /* Interpolate linearly between
                       current and next position */
                    if (i < N-1)
                        {
                        x1 = sorted_input[i].x;
                        x2 = sorted_input[i+1].x;
                        
                        quartile[j] = 
                            x1*( 1.5-k) +
                            x2*(-0.5+k);
                        /*
                          k=0.5 => quartile = x1,
                          k=1.0 => quartile = midpoint
                        */
                        }
                    else
                        {
                        quartile[j] = sorted_input[i].x;
                        }
                    }
                }
            }
        else
            {
            /* If there's no flux (sum = 0) then
               set quartiles to something that's not 
               completely insensible.
            */
            quartile[1] = cpl_matrix_get_mean(x_matrix);
            
            quartile[2] = quartile[1];
            quartile[0] = quartile[2] - 1.0;
            /* Then sigma_guess = 1.0 */
            }

        cpl_free(sorted_input);
        } /* If need to compute quartiles */
        
    /* x0_guess = median of distribution */
    x0_guess = (fit_pars & CPL_FIT_CENTROID) ? quartile[1] : *x0;
    
    /* sigma_guess = median of absolute residuals
     *
     *  (68% is within 1 sigma, and 50% is within 0.6744
     *  sigma,  => quartile3-quartile1 = 2*0.6744 sigma)
     */
    sigma_guess = (fit_pars & CPL_FIT_STDEV) ? 
        (quartile[2] - quartile[0]) / (2*0.6744) : *sigma;
    
    area_guess  = (fit_pars & CPL_FIT_AREA) ?
        (cpl_vector_get_max(y) - offset_guess) * sqrt(2*M_PI) * sigma_guess : *area;
    /* This formula makes sense only if the area is positive */
    
    /* Make sure that the area is a positive number */
    if ( area_guess <= 0)  area_guess = 1.0;
    if (sigma_guess <= 0) sigma_guess = 1.0;
    }
    
    /* Wrap parameters, fit, unwrap */
    {
    cpl_vector *a;
    int ia[5];            /* The last element
                 is ignored if
                 M = 4 */
    cpl_error_code ec;

    assert(M == 4 || M == 5);
    a = cpl_vector_new(M);

    if (a == NULL)
        {
        cpl_matrix_unwrap(x_matrix);
        cpl_ensure_code(
            CPL_FALSE,
            CPL_ERROR_ILLEGAL_OUTPUT);
        }

    cpl_vector_set(a, 0, x0_guess);
    cpl_vector_set(a, 1, sigma_guess);
    cpl_vector_set(a, 2, area_guess);
    cpl_vector_set(a, 3, offset_guess);
    
    ia[0] = fit_pars & CPL_FIT_CENTROID;
    ia[1] = fit_pars & CPL_FIT_STDEV;
    ia[2] = fit_pars & CPL_FIT_AREA;
    ia[3] = fit_pars & CPL_FIT_OFFSET;

    if (M == 5)
        {
        /* linear sky-term, first hold it constant,
         * then call LM-fitting a second time where
         * it is non-constant */
        if (fit_pars & CPL_FIT_OFFSET)
                    {
                        cpl_vector_set(a, 4, 0);
                    }
                else
                    {
                        cpl_vector_set(a, 4, *slope);
                    }

        ia[4] = 0;
        }
    
    ec = uves_fit(x_matrix, NULL,
            y, sigma_y, 
            a, ia, f, dfda,
            mse, red_chisq,
            covariance);
    
//    printf("Sky: e='%s'\n", cpl_error_get_message()); cpl_vector_dump(a, stdout);
    if (M == 5)
        {
        ia[4] = fit_pars & CPL_FIT_OFFSET;

        if (covariance != NULL)
            {
            uves_free_matrix(covariance);
            }

        ec = uves_fit(x_matrix, NULL,
                y, sigma_y, 
                a, ia, f, dfda,
                mse, red_chisq,
                covariance);
//    printf("Sky: e='%s'\n", cpl_error_get_message()); cpl_vector_dump(a, stdout);
        }

    cpl_matrix_unwrap(x_matrix);
    
    /* Check return status of fitting routine */
    if (ec == CPL_ERROR_NONE      ||
        ec == CPL_ERROR_SINGULAR_MATRIX)
        {
        /* The LM algorithm converged. The computation
         *  of the covariance matrix might have failed.
         */
        
        /* In principle, the LM algorithm might have converged
         * to a negative sigma (even if the guess value was
         * positive). Make sure that the returned sigma is positive
         * (by convention).
         */

        if (CPL_FIT_CENTROID) *x0     = cpl_vector_get(a, 0);
        if (CPL_FIT_STDEV   ) *sigma  = fabs(cpl_vector_get(a, 1));
        if (CPL_FIT_AREA    ) *area   = cpl_vector_get(a, 2);
        if (CPL_FIT_OFFSET  ) {
                    *offset = cpl_vector_get(a, 3);
                    if (M == 5) *slope = cpl_vector_get(a, 4);
                }
        }
    
    cpl_vector_delete(a);
    /* Min/max x                    */
    double xlo = cpl_vector_get_min(x);
    double xhi = cpl_vector_get_max(x);

    if (ec == CPL_ERROR_CONTINUE ||
        !(
        !irplib_isnan(*x0    ) && !irplib_isinf(*x0    ) &&
        !irplib_isnan(*sigma ) && !irplib_isinf(*sigma ) &&
        !irplib_isnan(*area  ) && !irplib_isinf(*area  ) &&
        !irplib_isnan(*offset) && !irplib_isinf(*offset) &&
        ((M != 5) || (!irplib_isnan(*slope ) && !irplib_isinf(*slope ))) &&
        xlo <= *x0 && *x0 <= xhi &&
        0 < *sigma && *sigma < (xhi - xlo + 1) &&
        (fit_back || 0 < *area)
        /* This extra check on the background level makes sense
           iff the input flux is assumed to be positive
           && *offset > - *area  */
        )
        )
        {
        /* The LM algorithm did not converge, or it converged to
         * a non-sensical result. Return the guess parameter values
         * in order to enable the caller to recover. */

        *x0         = x0_guess;
        *sigma      = sigma_guess;
        *area       = area_guess;
        *offset     = offset_guess;
                if (M == 5) *slope = 0;
        
        /* In this case the covariance matrix will not make sense
           (because the LM algorithm failed), so delete it */
        if (covariance != NULL && *covariance != NULL)
            {
            cpl_matrix_delete(*covariance);
            *covariance = NULL;
            }

        /* Return CPL_ERROR_CONTINUE if the fitting routine failed */
        cpl_ensure_code(
            CPL_FALSE,
            CPL_ERROR_CONTINUE);
        }
    }
    
    return CPL_ERROR_NONE;
}
/* @endcond */

#define DEBUG_LM 0   /* Set to non-zero to print info on the error msg level */
/*----------------------------------------------------------------------------*/
/**
   @brief   Fit a function to a set of data
   @param   x        N x D matrix of the positions to fit.
                     Each matrix row is a D-dimensional position.
   @param   sigma_x  Uncertainty (one sigma, gaussian errors assumed)
                     assosiated with @em x. Taking into account the 
             uncertainty of the independent variable is currently
             unsupported, and this parameter must therefore be set
             to NULL.
   @param   y        The N values to fit.
   @param   sigma_y  Vector of size N containing the uncertainties of
                     the y-values. If this parameter is NULL, constant
             uncertainties are assumed.
   @param   a        Vector containing M fit parameters. Must contain
                     a guess solution on input and contains the best
             fit parameters on output.
   @param   ia       Array of size M defining which fit parameters participate
                     in the fit (non-zero) and which fit parameters are held
             constant (zero). At least one element must be non-zero.
             Alternatively, pass NULL to fit all parameters.
   @param   f        Function that evaluates the fit function
                     at the position specified by the first argument (an array of
             size D) using the fit parameters specified by the second
             argument (an array of size M). The result must be output
             using the third parameter, and the function must return zero
             iff the evaluation succeded.
   @param   dfda     Function that evaluates the first order partial
                     derivatives of the fit function with respect to the fit
             parameters at the position specified by the first argument
             (an array of size D) using the parameters specified by the
             second argument (an array of size M). The result must
             be output using the third parameter (array of size M), and
             the function must return zero iff the evaluation succeded.
   @param mse        If non-NULL, the mean squared error of the best fit is
                     computed.
   @param red_chisq  If non-NULL, the reduced chi square of the best fit is
                     computed. This requires @em sigma_y to be specified.
   @param covariance If non-NULL, the formal covariance matrix of the best
                     fit parameters is computed (or NULL on error). On success
             the diagonal terms of the covariance matrix are guaranteed
             to be positive. However, terms that involve a constant
             parameter (as defined by the input array @em ia) are
             always set to zero. Computation of the covariacne matrix
             requires @em sigma_y to be specified.
            

   @return  CPL_ERROR_NONE iff OK.

   This function makes a minimum chi squared fit of the specified function
   to the specified data set using the Levenberg-Marquardt algorithm.

   Possible _cpl_error_code_ set in this function:
   - CPL_ERROR_NULL_INPUT if an input pointer other than @em sigma_x, @em
     sigma_y, @em mse, @em red_chisq or @em covariance is NULL.
   - CPL_ERROR_ILLEGAL_INPUT if an input matrix/vector is empty, if @em ia
     contains only zero values, if N <= M and @em red_chisq is non-NULL, 
     if any element of @em sigma_x or @em sigma_y is non-positive, or if
     evaluation of the fit function or its derivative failed.
   - CPL_ERROR_INCOMPATIBLE_INPUT if the dimensions of the input
     vectors/matrices do not match, or if chi square or covariance computation
     is requested and @em sigma_y is NULL.
   - CPL_ERROR_ILLEGAL_OUTPUT if memory allocation failed.
   - CPL_ERROR_CONTINUE if the Levenberg-Marquardt algorithm failed to converge.
   - CPL_ERROR_SINGULAR_MATRIX if the covariance matrix could not be computed.

*/
/*----------------------------------------------------------------------------*/
cpl_error_code
uves_fit(const cpl_matrix *x, const cpl_matrix *sigma_x,
     const cpl_vector *y, const cpl_vector *sigma_y,
     cpl_vector *a, const int ia[],
     int    (*f)(const double x[], const double a[], double *result),
     int (*dfda)(const double x[], const double a[], double result[]),
     double *mse,
     double *red_chisq,
     cpl_matrix **covariance)
{
    const double *x_data     = NULL; /* Pointer to input data                  */
    const double *y_data     = NULL; /* Pointer to input data                  */
    const double *sigma_data = NULL; /* Pointer to input data                  */
    int N    = 0;                    /* Number of data points                  */
    int D    = 0;                    /* Dimension of x-points                  */
    int M    = 0;                    /* Number of fit parameters               */
    int Mfit = 0;                    /* Number of non-constant fit
                        parameters                             */

    double lambda    = 0.0;          /* Lambda in L-M algorithm                */
    double MAXLAMBDA = 10e40;        /* Parameter to control the graceful exit
                    if steepest descent unexpectedly fails */
    double chi_sq    = 0.0;          /* Current  chi^2                         */
    int count        = 0;            /* Number of successive small improvements
                    in chi^2 */
    int iterations   = 0;
   
    cpl_matrix *alpha  = NULL;       /* The MxM ~curvature matrix used in L-M  */
    cpl_matrix *beta   = NULL;       /* Mx1 matrix = -.5 grad(chi^2)           */
    double *a_data     = NULL;       /* Parameters, a                          */
    double *a_da       = NULL;       /* Candidate position a+da                */
    double *part       = NULL;       /* The partial derivatives df/da          */
    int *ia_local      = NULL;       /* non-NULL version of ia                 */
   
    /* If covariance computation is requested, then either
     * return the covariance matrix or return NULL.
     */
    if (covariance != NULL) *covariance = NULL;

    /* Validate input */
    cpl_ensure_code(x       != NULL, CPL_ERROR_NULL_INPUT);
    cpl_ensure_code(sigma_x == NULL, CPL_ERROR_UNSUPPORTED_MODE);
    cpl_ensure_code(y       != NULL, CPL_ERROR_NULL_INPUT);
    cpl_ensure_code(a       != NULL, CPL_ERROR_NULL_INPUT);
    /* ia may be NULL */
    cpl_ensure_code(f       != NULL, CPL_ERROR_NULL_INPUT);
    cpl_ensure_code(dfda    != NULL, CPL_ERROR_NULL_INPUT);

    /* Chi^2 and covariance computations require sigmas to be known */
    cpl_ensure_code( sigma_y != NULL || (red_chisq == NULL && covariance == NULL),
             CPL_ERROR_INCOMPATIBLE_INPUT);

    D = cpl_matrix_get_ncol(x);
    N = cpl_matrix_get_nrow(x);
    M = cpl_vector_get_size(a);
    cpl_ensure_code(N > 0 && D > 0 && M > 0, CPL_ERROR_ILLEGAL_INPUT);

    cpl_ensure_code( cpl_vector_get_size(y) == N,
             CPL_ERROR_INCOMPATIBLE_INPUT);

    x_data = cpl_matrix_get_data_const(x);
    y_data = cpl_vector_get_data_const(y);
    a_data = cpl_vector_get_data(a);

    if (sigma_y != NULL)
    {
        cpl_ensure_code( cpl_vector_get_size(sigma_y) == N,
                 CPL_ERROR_INCOMPATIBLE_INPUT);
        /* Sigmas must be positive */
        cpl_ensure_code( cpl_vector_get_min (sigma_y) > 0,
                 CPL_ERROR_ILLEGAL_INPUT);
        sigma_data = cpl_vector_get_data_const(sigma_y);
    }

    ia_local = cpl_malloc(M * sizeof(int));
    cpl_ensure_code(ia_local != NULL, CPL_ERROR_ILLEGAL_OUTPUT);

    /* Count non-constant fit parameters, copy ia */
    if (ia != NULL)
    {
        int i;

        Mfit = 0;
        for (i = 0; i < M; i++)
        {
            ia_local[i] = ia[i];

            if (ia[i] != 0) 
            {
                Mfit += 1;
            }
        }
        
        if (! (Mfit > 0))
        {
            cpl_free(ia_local);
            cpl_ensure_code( CPL_FALSE,
                     CPL_ERROR_ILLEGAL_INPUT);
        }
    }
    else
    {
        /* All parameters participate */
        int i;
        
        Mfit = M;
        
        for (i = 0; i < M; i++)
        {
            ia_local[i] = 1;
        }
    }

    /* To compute reduced chi^2, we need N > Mfit */
    if (! ( red_chisq == NULL || N > Mfit ) )
    {
        cpl_free(ia_local);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_INPUT);
    }

    /* Create alpha, beta, a_da, part  work space */
    alpha = cpl_matrix_new(Mfit, Mfit);
    if (alpha == NULL)
    {
        cpl_free(ia_local);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_OUTPUT);
    }
   
    beta = cpl_matrix_new(Mfit, 1);
    if (beta == NULL)
    {
        cpl_free(ia_local);
        cpl_matrix_delete(alpha);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_OUTPUT);
    }

    a_da = cpl_malloc(M * sizeof(double));
    if (a_da == NULL)
    {
        cpl_free(ia_local);
        cpl_matrix_delete(alpha);
        cpl_matrix_delete(beta);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_OUTPUT);
    }

    part = cpl_malloc(M * sizeof(double));
    if (part == NULL)
    {
        cpl_free(ia_local);
        cpl_matrix_delete(alpha);
        cpl_matrix_delete(beta);
        cpl_free(a_da);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_ILLEGAL_OUTPUT);
    }

    /* Initialize loop variables */
    lambda = 0.001;
    count = 0;
    iterations = 0;
    if( (chi_sq = get_chisq(N, D, f, a_data, x_data, y_data, sigma_data)) < 0)
    {
        cpl_free(ia_local);
        cpl_matrix_delete(alpha);
        cpl_matrix_delete(beta);
        cpl_free(a_da);
        cpl_free(part);
        cpl_ensure_code(
        CPL_FALSE,
        cpl_error_get_code());
    }

#if DEBUG_LM    
     uves_msg_error("Initial chi^2 = %f", chi_sq); 
     {int i;
     for (i = 0; i < M; i++) 
     {
         uves_msg_error("Initial a[%d] = %e", i, a_data[i]); 
     }
     }
#endif

    /* Iterate until chi^2 didn't improve substantially many (say, 5)
       times in a row */
    while (count < 5 && lambda < MAXLAMBDA && iterations < UVES_FIT_MAXITER)
    {
        /* In each iteration lambda increases, or chi^2 decreases or
           count increases. Because chi^2 is bounded from below
           (and lambda and count from above), the loop will terminate */

        double chi_sq_candidate = 0.0;
        int returncode = 0;

        /* Get candidate position in parameter space = a+da,
         * where  alpha * da = beta .
         * Increase lambda until alpha is non-singular
         */
       
        while( (returncode = get_candidate(a_data, ia_local,
                           M, N, D,
                           lambda, f, dfda,
                           x_data, y_data, sigma_data,
                           part, alpha, beta, a_da)
               ) != 0
           && cpl_error_get_code() == CPL_ERROR_SINGULAR_MATRIX
           && lambda < MAXLAMBDA)
        {
#if DEBUG_LM    
            uves_msg_error("Singular matrix. lambda = %e", lambda);
#endif
            cpl_error_reset();
            lambda *= 9.0;
        }
       
        /* Set error if lambda diverged */
        if ( !( lambda < MAXLAMBDA ) )
        {
            cpl_free(ia_local);
            cpl_matrix_delete(alpha);
            cpl_matrix_delete(beta);
            cpl_free(a_da);
            cpl_free(part);
            cpl_ensure_code(
            CPL_FALSE,
            CPL_ERROR_CONTINUE);
        }
       
        if (returncode != 0)
        {
            cpl_free(ia_local);
            cpl_matrix_delete(alpha);
            cpl_matrix_delete(beta);
            cpl_free(a_da);
            cpl_free(part);
            cpl_ensure_code(
            CPL_FALSE,
            cpl_error_get_code());
        }

        /* Get chi^2(a+da) */
        if ((chi_sq_candidate = get_chisq(N, D, f, a_da,
                          x_data, y_data, sigma_data)) < 0)
        {
            cpl_free(ia_local);
            cpl_matrix_delete(alpha);
            cpl_matrix_delete(beta);
            cpl_free(a_da);
            cpl_free(part);
            cpl_ensure_code(
            CPL_FALSE,
            cpl_error_get_code());
        }

        if (chi_sq_candidate > 1.000001 * chi_sq)
        {
            /* Move towards steepest descent */
#if DEBUG_LM
            uves_msg_error("Chi^2 = %f  Candidate = %f  Lambda = %e",
               chi_sq, chi_sq_candidate, lambda); 
#endif            
            lambda *= 9.0;
        }
        else
        {
#if DEBUG_LM
            uves_msg_error("Chi^2 = %f  Candidate = %f* Lambda = %e count = %d",
               chi_sq, chi_sq_candidate, lambda, count);
#endif
       
            /* Move towards Newton's algorithm */
            lambda /= 10.0;

            /* Count the number of successive improvements in chi^2 of
               less than 0.01 relative */
            if ( chi_sq == 0 ||
             (chi_sq - chi_sq_candidate)/chi_sq < .01)
            {
                count += 1;
            }
            else
            {
                /* Chi^2 improved by a significant amount,
                   reset counter */
                count = 0;
            }

            if (chi_sq_candidate < chi_sq)
            {
                /* chi^2 improved, update a and chi^2 */
                int i;
                for (i = 0; i < M; i++) 
                {
                    a_data[i] = a_da[i];
#if DEBUG_LM
                    uves_msg_error("-> a[%d] = %e", i, a_da[i]); 
#endif
                }
                chi_sq = chi_sq_candidate;
            }
        }
        iterations++;
    }

    /* Set error if we didn't converge */
    if ( !( lambda < MAXLAMBDA && iterations < UVES_FIT_MAXITER ) )
    {
#if DEBUG_LM
        uves_msg_error("Failed to converge, lambda=%f iterations=%d",
               lambda, iterations);
#endif
        cpl_free(ia_local);
        cpl_matrix_delete(alpha);
        cpl_matrix_delete(beta);
        cpl_free(a_da);
        cpl_free(part);
        cpl_ensure_code(
        CPL_FALSE,
        CPL_ERROR_CONTINUE);
    }

    /* Compute mse if requested */
    if (mse != NULL)
    {
        int i;

        *mse = 0.0;
       
        for(i = 0; i < N; i++)
        {
            double fx_i = 0.0;
            double residual = 0.0;
           
            /* Evaluate f(x_i) at the best fit parameters */
            if( f(&(x_data[i*D]),
              a_data,
              &fx_i) != 0)
            {
                cpl_free(ia_local);
                cpl_matrix_delete(alpha);
                cpl_matrix_delete(beta);
                cpl_free(a_da);
                cpl_free(part);
                cpl_ensure_code(
                CPL_FALSE,
                CPL_ERROR_ILLEGAL_INPUT);
            }

            residual = y_data[i] - fx_i;
            *mse += residual * residual;
        }
        *mse /= N;
    }

    /* Compute reduced chi^2 if requested */
    if (red_chisq != NULL)
    {
        /* We already know the optimal chi^2 (and that N > Mfit)*/
        *red_chisq = chi_sq / (N-Mfit);
    }

    /* Compute covariance matrix if requested
     * cov = alpha(lambda=0)^-1              
     */
    if (covariance != NULL)
    {
        cpl_matrix *cov;

        if( get_candidate(a_data, ia_local, 
                  M, N, D, 0.0, f, dfda, 
                  x_data, y_data, sigma_data,
                  part, alpha, beta, a_da)
        != 0)
        {
            cpl_free(ia_local);
            cpl_matrix_delete(alpha);
            cpl_matrix_delete(beta);
            cpl_free(a_da);
            cpl_free(part);
            cpl_ensure_code(
            CPL_FALSE,
            cpl_error_get_code());
        }
       
        cov = cpl_matrix_invert_create(alpha);
        if (cov == NULL)
        {
            cpl_free(ia_local);
            cpl_matrix_delete(alpha);
            cpl_matrix_delete(beta);
            cpl_free(a_da);
            cpl_free(part);
            cpl_ensure_code(
            CPL_FALSE,
            cpl_error_get_code());
        }
       
        /* Make sure that variances are positive */
        {
        int i;
        for (i = 0; i < Mfit; i++)
            {
            if ( !(cpl_matrix_get(cov, i, i) > 0) )
                {
                cpl_free(ia_local);
                cpl_matrix_delete(alpha);
                cpl_matrix_delete(beta);
                cpl_free(a_da);
                cpl_free(part);
                cpl_matrix_delete(cov);
                *covariance = NULL;
                cpl_ensure_code(
                    CPL_FALSE,
                    CPL_ERROR_SINGULAR_MATRIX);
                }
            }
        }

        /* Expand covariance matrix from Mfit x Mfit
           to M x M. Set rows/columns corresponding to fixed
           parameters to zero */

        *covariance = cpl_matrix_new(M, M);
        if (*covariance == NULL)
        {
            cpl_free(ia_local);
            cpl_matrix_delete(alpha);
            cpl_matrix_delete(beta);
            cpl_free(a_da);
            cpl_free(part);
            cpl_matrix_delete(cov);
            cpl_ensure_code(
            CPL_FALSE,
            CPL_ERROR_ILLEGAL_OUTPUT);
        }

        {
        int j, jmfit;

        for (j = 0, jmfit = 0; j < M; j++)
            if (ia_local[j] != 0)
            {
                int i, imfit;

                for (i = 0, imfit = 0; i < M; i++)
                if (ia_local[i] != 0)
                    {
                    cpl_matrix_set(*covariance, i, j,
                               cpl_matrix_get(
                               cov, imfit, jmfit));
                    imfit += 1;
                    }
                
                assert( imfit == Mfit );

                jmfit += 1;
            }
        
        assert( jmfit == Mfit );
        }

        cpl_matrix_delete(cov);
    }
    
    cpl_free(ia_local);
    cpl_matrix_delete(alpha);
    cpl_matrix_delete(beta);
    cpl_free(a_da);
    cpl_free(part);
   
    return CPL_ERROR_NONE;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Cast an image
   @param    image           Image to cast
   @param    to_type         The image is casted to this type
   @return   CPL_ERROR_NONE iff OK

   This function wraps @c cpl_image_cast().

*/
/*----------------------------------------------------------------------------*/
cpl_error_code
uves_cast_image(cpl_image **image, cpl_type to_type)
{
    cpl_image *temp = NULL;
    
    assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");

    temp = *image;
    
    check( *image = cpl_image_cast(temp, to_type), "Couldn't convert image to %s", 
       uves_tostring_cpl_type(to_type));
    
  cleanup:
    uves_free_image(&temp);
    return cpl_error_get_code();    
}


/*-----------------------------------------------------------------------------*/
/**
   @brief   Crop an image
   @param   image   Image to crop
   @param   x1      Lower left (inclusive)
   @param   y_1     Lower left (inclusive)
   @param   x2      Upper right (inclusive)
   @param   y2      Upper right (inclusive)
   @return  CPL_ERROR_NONE iff OK

   This function wraps @c cpl_image_extract().

   The sub-window (@em x1, @em y1) - (@em x2, @em y2) is extracted.
   
*/
/*-----------------------------------------------------------------------------*/
cpl_error_code
uves_crop_image(cpl_image **image, int x1, int y_1, int x2, int y2)
{
    cpl_image *temp = NULL;
    
    assure( image != NULL, CPL_ERROR_NULL_INPUT, "Null image");

    temp = *image;
    
    check( *image = cpl_image_extract(temp, x1, y_1, x2, y2), 
       "Could not extract image");
    
  cleanup:
    uves_free_image(&temp);
    return cpl_error_get_code();    
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Read a property value from a property list
   @param    plist       Propertylist to read
   @param    keyword     Name of property to read
   @param    keywordtype Type of keyword
   @param    result      The value read

   @return   CPL_ERROR_NONE iff OK

   This function wraps  @c uves_propertylist_get_int(), 
   @c uves_propertylist_get_bool(), @c uves_propertylist_get_double()
   and @c uves_propertylist_get_string(). It checks existence and type of the
   requested keyword before reading, and gives informative error messages
   if the property could not be read.

   @note The result is written to the variable pointed to by the parameter
   @em result. Because this is a void pointer, it is the responsibility of
   the caller to make sure that the type of this pointer variable corresponds
   to the requested @em keywordtype. E.g. if @em keywordtype is CPL_TYPE_INT,
   then @em result must be an int pointer (int *). If @em keywordtype is
   CPL_TYPE_STRING, then @em result must be a char **, and so on.

*/
/*----------------------------------------------------------------------------*/
cpl_error_code
uves_get_property_value(const uves_propertylist *plist, const char *keyword, 
            cpl_type keywordtype, void *result)
{
    cpl_type t;
    
    /* Check input */
    assure( plist != NULL  , CPL_ERROR_NULL_INPUT, "Null property list");
    assure( keyword != NULL, CPL_ERROR_NULL_INPUT, "Null keyword");
    
    /* Check for existence... */
    assure( uves_propertylist_contains(plist, keyword), CPL_ERROR_DATA_NOT_FOUND,
        "Keyword %s does not exist", keyword );
    /* ...and type of keyword */
    check( t = uves_propertylist_get_type(plist, keyword) ,
       "Could not read type of keyword '%s'", keyword );
    assure(t == keywordtype , CPL_ERROR_TYPE_MISMATCH, 
       "Keyword '%s' has wrong type (%s). %s expected",
       keyword, uves_tostring_cpl_type(t), uves_tostring_cpl_type(keywordtype));
    
    /* Read the keyword */
    switch (keywordtype)
    {
    case CPL_TYPE_INT   : 
        check( *((      int    *)result) =
        uves_propertylist_get_int(plist, keyword),
        "Could not get (integer) value of %s", keyword );
        break;
    case CPL_TYPE_BOOL  : 
        check( *((      bool   *)result) =
           uves_propertylist_get_bool(plist, keyword), 
           "Could not get (boolean) value of %s", keyword ); 
        break;
    case CPL_TYPE_DOUBLE:
        check( *((      double *)result) = 
           uves_propertylist_get_double(plist, keyword), 
           "Could not get (double) value of %s" , keyword );
        break;
    case CPL_TYPE_STRING: 
        check( *((const char * *)result) = 
           uves_propertylist_get_string(plist, keyword), 
           "Could not get (string) value of %s" , keyword ); 
        break;
    default:
        assure( false, CPL_ERROR_INVALID_TYPE, "Unknown type");
        break;
    }
    
  cleanup:
    return cpl_error_get_code();
}


/*----------------------------------------------------------------------------*/
/**
   @brief    Find a frame in a frame set
   @param    frames     Frame set to search
   @param    wanted     List of requested tags
   @param    N          Size of @em wanted
   @param    found      (output) index of identified tag
   @param    frame      (output) if non-NULL, this will point to the identified
                        frame (or to NULL on error)

   @return   Filename of frame, or NULL if no matching frame was found

   This function wraps  @c cpl_frameset_find(), but searches for a frame 
   matching one of several tags.

   The parameter @em wanted is an array of tags to search for. Upon termination
   the parameter @em found will contain the index of the identified tag. The
   filename associated with the identified frame is returned. If no matching
   frame could be found, an error is set, NULL is returned, and @em found is undefined.

   @note The function will try to read the first @em N elements of the array
   @em wanted, so it is the responsibility of the caller that @em N is really
   the size of the array

*/
/*----------------------------------------------------------------------------*/

const char *
uves_find_frame(const cpl_frameset *frames, const char **wanted, int N, int *found,
        const cpl_frame **frame)
{
    const char *filename = NULL;  /* Return NULL in case of error */
    int i;
    const cpl_frame *f = NULL;
    
    /* Return well-defined pointers in case of error */
    *found = 0;
    if (frame != NULL)
    {
        *frame = NULL;
    }

    for (i = 0; i < N; i++)
    {
        check( f = cpl_frameset_find_const(frames, wanted[i]), 
           "Could not search through frame set");
        if (f != NULL) 
        {
            check( filename = cpl_frame_get_filename(f), 
               "Could not read frame filename");
            *found = i;
            if (frame != NULL)
            {
                *frame = f;
            }
            /* break */
            i = N;
        }
    }
    
    /* Set an error if a matching frame wasn't found */
    assure(filename != NULL, CPL_ERROR_DATA_NOT_FOUND, "No such frame in frame set");

  cleanup:
    return filename;
}


/*----------------------------------------------------------------------------*/
/**
   @brief    Get number of extensions of a FITS file
   @param    filename      The FITS file
   @return   The number of extensions in the specified file, or undefined on error

*/
/*----------------------------------------------------------------------------*/
cpl_size
uves_get_nextensions(const char *filename)
{
    cpl_size result = 0;
    cpl_frame *f = NULL;
    
    /* CPL only supports reading the number of extensions of a FITS
       file if this is in a frame, so create a frame for the specified file */
      
    check(( f = cpl_frame_new(),
        cpl_frame_set_filename(f, filename)),
      "Could not create frame");

    check( result = cpl_frame_get_nextensions(f),
       "Error reading number of extensions of file '%s'", filename);
  cleanup:
    cpl_frame_delete(f);
    return result;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Read a parameter value from a parameter list
   @param    parameters     A parameter list
   @param    context        Context of parameter within calling recipe
   @param    recipe_id      Name of calling recipe
   @param    name           Full name of parameter to read is @em context . @em name
   @param    type           Type of parameter to read
   @param    value          (output) Value of the parameter read
   @return   CPL_ERROR_NONE iff OK

   This function does error checking (existence, type) before trying
   to read the parameter.

   @note The @em value parameter is passed as a void pointer. Therefore, it 
   is the responsibility of the caller to make sure that the type of this
   pointer corresponds to the requested parameter type. E.g. if @em type 
   is CPL_TYPE_INT, then @em value must be an int pointer (int *). If 
   @em keywordtype is CPL_TYPE_STRING, then @em result must be a char **,
   and so on.
*/
/*----------------------------------------------------------------------------*/
cpl_error_code
uves_get_parameter(const cpl_parameterlist *parameters, const char *context,
           const char *recipe_id, const char *name, cpl_type type, 
           void *value)
{
    char *fullname = NULL;
    const cpl_parameter *p = NULL;
    
    passure( parameters != NULL, " ");
    /* 'context' may be NULL */
    passure( recipe_id != NULL, " ");
    passure( name != NULL, " ");
    passure( value != NULL, " ");

    if (context != NULL)
    {
        check( fullname = uves_sprintf("%s.%s.%s", context, recipe_id, name),
           "Error getting full parameter name");
    }
    else
    {
        check( fullname = uves_sprintf("%s.%s", recipe_id, name),
           "Error getting full parameter name");
    }

    /* Const cast */
    check( p = cpl_parameterlist_find_const(parameters, fullname), 
       "Error searching for parameter '%s'", fullname);
    assure( p != NULL, CPL_ERROR_DATA_NOT_FOUND, 
        "No parameter '%s' in parameter list", fullname);
    
    /* Check that parameter has the correct type */
    {
    cpl_type partype;
    check(  partype = cpl_parameter_get_type(p), 
        "Could not read type of parameter '%s'", fullname);
    assure( partype == type, CPL_ERROR_TYPE_MISMATCH,
        "Parameter '%s' has type %s. Expected type was %s", 
        fullname,
        uves_tostring_cpl_type(partype), uves_tostring_cpl_type(type));    
    }

    /* Read the parameter */
    switch(type)
    {
    case CPL_TYPE_INT:
        check( *(int *         )value = 
           cpl_parameter_get_int   (p),
           "Could not read integer parameter '%s'", fullname);  
        break;
    case CPL_TYPE_BOOL:
        check( *(bool *        )value = 
           cpl_parameter_get_bool  (p),
           "Could not read boolean parameter '%s'", fullname);  
        break;
    case CPL_TYPE_DOUBLE:    
        check( *(double *      )value = 
           cpl_parameter_get_double(p),
           "Could not read double parameter '%s'" , fullname );  
        break;
    case CPL_TYPE_STRING:
        check( *(const char ** )value = 
           cpl_parameter_get_string(p),
           "Could not read string parameter '%s'" , fullname );  
        break;
    default:
        assure(false, CPL_ERROR_UNSUPPORTED_MODE,
           "Don't know how to read parameter '%s' of type %s",
           fullname, uves_tostring_cpl_type(type));
        break;
    }
    
  cleanup:
    cpl_free(fullname); fullname = NULL;
    return cpl_error_get_code();
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Set a parameter
   @param    parameters     A parameter list
   @param    context        Context of parameter
   @param    name           Full name of parameter to read is @em context . @em name
   @param    type           Type of parameter
   @param    value          The value to set
   @return   CPL_ERROR_NONE iff OK
   
   This function does error checking (existence, type) before setting the parameter.

   @note The @em value parameter is passed as a void pointer. This parameter 
   must be the address of a variable (of appropriate type) containing the 
   value to set.
*/
/*----------------------------------------------------------------------------*/
cpl_error_code
uves_set_parameter(cpl_parameterlist *parameters, 
           const char *context,
           const char *name, cpl_type type, void *value)
{
    char *fullname = NULL;
    cpl_parameter *p = NULL;
    
    check( fullname = uves_sprintf("%s.%s", context, name),
       "Error getting full parameter name");

    /* Const cast */
    check( p = cpl_parameterlist_find(parameters, fullname), 
       "Error searching for parameter '%s'", fullname);
    assure( p != NULL, CPL_ERROR_DATA_NOT_FOUND, 
        "No parameter '%s' in parameter list", fullname);

    /* Check that parameter has the correct type */
    {
    cpl_type partype;
    check(  partype = cpl_parameter_get_type(p), 
        "Could not read type of parameter '%s'", fullname);
    assure( partype == type, CPL_ERROR_TYPE_MISMATCH,
        "Parameter '%s' has type %s. Expected type was %s", 
        fullname, uves_tostring_cpl_type(partype),
        uves_tostring_cpl_type(type));    
    }

    /* Set the parameter */
    switch(type)
    {
    case CPL_TYPE_INT:
        check( cpl_parameter_set_int   (p, *((int *)    value)), 
           "Could not set integer parameter '%s'", fullname);  
        break;
    case CPL_TYPE_BOOL:    
        check( cpl_parameter_set_bool  (p, *((bool *)    value)), 
           "Could not set boolean parameter '%s'", fullname);
        break;
    case CPL_TYPE_DOUBLE:
        check( cpl_parameter_set_double(p, *((double *) value)), 
           "Could not set double parameter '%s'" , fullname);
        break;
    case CPL_TYPE_STRING:
        check( cpl_parameter_set_string(p, *((char **)  value)),
           "Could not set string parameter '%s'" , fullname);  
        break;
    default:
        assure(false, CPL_ERROR_UNSUPPORTED_MODE,
           "Don't know how to set parameter of type %s", 
           uves_tostring_cpl_type(type));
        break;
    }
    
  cleanup:
    cpl_free(fullname); fullname = NULL;
    return cpl_error_get_code();
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Change the default value of a parameter
   @param    parameters     A parameter list
   @param    context        context of parameter, or NULL
   @param    parname        Name of parameter to set
   @param    type           Type of parameter
   @param    value          The new default value to set
   @return   CPL_ERROR_NONE iff OK

   This function provides a minmal interface for changing a parameter's
   default value, given a parameter list.

*/
/*----------------------------------------------------------------------------*/

cpl_error_code
uves_set_parameter_default(cpl_parameterlist *parameters, const char *context,
               const char *parname, 
               cpl_type type, void *value)
{
    const char *full_name = NULL;
    cpl_parameter *p = NULL;
    cpl_type partype;

    if (context != NULL)
    {
        full_name = uves_sprintf("%s.%s", context, parname);
    }
    else
    {
        full_name = uves_sprintf("%s", parname);
    }

    if (full_name == NULL)
    {
        return CPL_ERROR_ILLEGAL_OUTPUT;
    }


    if ( (p = cpl_parameterlist_find(parameters, full_name)) == NULL)
    {
        uves_msg_error("Missing parameter: '%s'", full_name);

        uves_free_string_const(&full_name);
        if (cpl_error_get_code() != CPL_ERROR_NONE) return cpl_error_get_code();
        else return CPL_ERROR_DATA_NOT_FOUND;
    }
    
    partype = cpl_parameter_get_type(p);
    
    if (partype != type)
    {
        uves_msg_error("Parameter '%s' has type %s. Expected type was %s", 
              full_name, uves_tostring_cpl_type(partype), 
              uves_tostring_cpl_type(type));

        uves_free_string_const(&full_name);
        return CPL_ERROR_TYPE_MISMATCH;
    }

    switch(type)
    {
    case CPL_TYPE_INT:
        cpl_parameter_set_default_int   (p, *((int *)    value)); 
        break;
    case CPL_TYPE_BOOL:
        cpl_parameter_set_default_bool  (p, *((bool *)   value)); 
        break;
    case CPL_TYPE_DOUBLE:
        cpl_parameter_set_default_double(p, *((double *) value)); 
        break;
    case CPL_TYPE_STRING:
        cpl_parameter_set_default_string(p, *((char **)  value)); 
        break;
    default:
        uves_msg_error("Unknown type: %s", uves_tostring_cpl_type(type));

        uves_free_string_const(&full_name);
        return CPL_ERROR_INVALID_TYPE;
    }

    if (cpl_error_get_code() != CPL_ERROR_NONE)
    {
        uves_msg_error("Error changing value of parameter '%s'", 
               full_name);

        uves_free_string_const(&full_name);
        return cpl_error_get_code();
    }


    uves_free_string_const(&full_name);
    return CPL_ERROR_NONE;
}


/*----------------------------------------------------------------------------*/
/**
   @brief    Change values lower than median to fractin of median
   @param    t        Table
   @param    column   Column name, must have type CPL_TYPE_DOUBLE
   @param    fraction fraction of median

   The values in the specified table @em column that are lower than the given
   fraction of the median of all column values are changed to that value
   (thresholding).

*/
/*----------------------------------------------------------------------------*/
void 
uves_raise_to_median_frac(cpl_table *t, const char *column, double fraction)
{
    int i = 0;
    double threshold;

    assure_nomsg(t != NULL, CPL_ERROR_NULL_INPUT);
    assure(cpl_table_has_column(t, column), CPL_ERROR_DATA_NOT_FOUND,
       "No such column: %s", column);

    assure(cpl_table_get_column_type(t, column) == CPL_TYPE_DOUBLE,
       CPL_ERROR_UNSUPPORTED_MODE, "Column %s has type %s. %s expected",
       column,
       uves_tostring_cpl_type(cpl_table_get_column_type(t, column)),
       uves_tostring_cpl_type(CPL_TYPE_DOUBLE));


    threshold = fraction * cpl_table_get_column_median(t, column);
    for (i = 0; i < cpl_table_get_nrow(t); i++)
    {
        if (cpl_table_get_double(t, column, i, NULL) < threshold)
        {
            cpl_table_set_double(t, column, i, threshold);
        }
    }

  cleanup:
    return;
}


/*----------------------------------------------------------------------------*/
/**
   @brief    Select table rows
   @param    t        Table
   @param    column   Column name
   @param    operator Logical operator
   @param    value    Value used for comparison
   @return   Number of selected rows

   A row is selected if and only if the value in @em column is in the
   relation @em operator to the specified @em value. The specified column must
   have type CPL_TYPE_DOUBLE, CPL_TYPE_FLOAT or CPL_TYPE_INT. If integer,
   the integer nearest to @em value is used for the comparison.

   Also see @c cpl_table_and_selected_<type>().

*/
/*----------------------------------------------------------------------------*/

int
uves_select_table_rows(cpl_table *t,  const char *column, 
               cpl_table_select_operator operator, double value)
{
    int result = 0;
    cpl_type type;
    
    assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT, 
        "No such column: %s", column);

    type = cpl_table_get_column_type(t, column);

    assure( type == CPL_TYPE_DOUBLE || type == CPL_TYPE_FLOAT ||
        type == CPL_TYPE_INT, CPL_ERROR_INVALID_TYPE,
        "Column '%s' must be double or int. %s found", column, 
        uves_tostring_cpl_type(type));

    check( cpl_table_select_all(t), "Error selecting rows");
    
    if      (type == CPL_TYPE_DOUBLE)
    {
        result = cpl_table_and_selected_double(t, column, operator, value);
    }
    else if (type == CPL_TYPE_FLOAT)
    {
        result = cpl_table_and_selected_float(t, column, operator, value);
    }
    else if (type == CPL_TYPE_INT)
    {
        result = cpl_table_and_selected_int(t, column, operator, 
                                                uves_round_double(value));
    }
    else { /*impossible*/ passure(false, ""); }
    
  cleanup:
    return result;

}

/*----------------------------------------------------------------------------*/
/**
   @brief    Extract table rows
   @param    t        Table
   @param    column   Column name
   @param    operator Logical operator
   @param    value    Value used for comparison
   @return   Number of rows left

   This function is just like @c uves_erase_table_rows, except
   the logic is reversed (only the rows matching the criterion
   are kept).

*/
/*----------------------------------------------------------------------------*/
int
uves_extract_table_rows_local(cpl_table *t, const char *column,
                  cpl_table_select_operator operator, double value)
{
    int result = 0;
    
    assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT, 
        "No such column: %s", column);

    check( result = uves_select_table_rows(t, column, operator, value),
       "Error selecting rows");

    check( cpl_table_not_selected(t), "Error selecting rows");

    check( cpl_table_erase_selected(t), "Error deleting rows");

  cleanup:
    return result;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Erase table rows
   @param    t        Table
   @param    column   Column name
   @param    operator Logical operator
   @param    value    Value used for comparison
   @return   Number of erased rows

   A table row is erased if and only if the value in @em column is in the
   relation @em operator to the specified @em value. The specified column must
   have type CPL_TYPE_DOUBLE or CPL_TYPE_INT. If integer, the integer nearest
   to @em value is used for the comparison.

   The table selection flags are reset.

   Also see @c cpl_table_and_selected_<type>().

*/
/*----------------------------------------------------------------------------*/
int
uves_erase_table_rows(cpl_table *t, const char *column,
              cpl_table_select_operator operator, double value)
{
    int result = 0;
    
    assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT, 
        "No such column: %s", column);

    check( result = uves_select_table_rows(t, column, operator, value),
       "Error selecting rows");

    check( cpl_table_erase_selected(t), "Error deleting rows");

  cleanup:
    return result;
}





/*----------------------------------------------------------------------------*/
/**
   @brief    append property to propertylist
   @param    plist    propertylist
   @param    p        property to append

   The comment is not copied

*/
/*----------------------------------------------------------------------------*/

void uves_propertylist_append_property(uves_propertylist *plist, const cpl_property *p)
{
    switch(cpl_property_get_type(p)) {
    case CPL_TYPE_CHAR:
        uves_propertylist_append_char(plist, cpl_property_get_name(p), cpl_property_get_char(p));
        break;
    case CPL_TYPE_BOOL:
        uves_propertylist_append_bool(plist, cpl_property_get_name(p), cpl_property_get_bool(p));
        break;
    case CPL_TYPE_INT:
        uves_propertylist_append_int(plist, cpl_property_get_name(p), cpl_property_get_int(p));
        break;
    case CPL_TYPE_LONG:
        uves_propertylist_append_long(plist, cpl_property_get_name(p), cpl_property_get_long(p));
        break;
    case CPL_TYPE_FLOAT:
        uves_propertylist_append_float(plist, cpl_property_get_name(p), cpl_property_get_float(p));
        break;
    case CPL_TYPE_DOUBLE:
        uves_propertylist_append_double(plist, cpl_property_get_name(p), cpl_property_get_double(p));
        break;
    case CPL_TYPE_STRING:
        uves_propertylist_append_string(plist, cpl_property_get_name(p), cpl_property_get_string(p));
        break;
    default:
        assure( false, CPL_ERROR_UNSUPPORTED_MODE,
                "Type is %s", uves_tostring_cpl_type(cpl_property_get_type(p)));
        break;
    }
  cleanup:
    return;
}



/*----------------------------------------------------------------------------*/
/**
   @brief    Workaround for broken CPL function
   @param    t        Table
   @param    column   Column name
*/
/*----------------------------------------------------------------------------*/
int
uves_table_and_selected_invalid(cpl_table *t, const char *column)
{
    if (cpl_table_get_column_type(t, column) != CPL_TYPE_STRING)
        {
            return cpl_table_and_selected_invalid(t, column);
        }
    else
        {
            int i = 0;
            for (i = 0; i < cpl_table_get_nrow(t); i++)
                {
                    if (cpl_table_is_selected(t, i))
                        {
                            if (cpl_table_is_valid(t, column, i))
                                {
                                    cpl_table_unselect_row(t, i);
                                }
                            /* else keep it selected */
                        }
                    /* else unselected, don't change */
                }

            return cpl_table_count_selected(t);
        }
}

    
/*----------------------------------------------------------------------------*/
/**
   @brief    Erase invalid table rows
   @param    t        Table
   @param    column   Column name or NULL
   @return   Number of erased rows

   The function erases rows where the value of @em column is invalid.

   If @em column is NULL, rows with invalid value in any column are erased.

   No column is ever removed, even if it contains only invalid elements (which 
   is why this function is different from @c cpl_table_erase_invalid()).

   The table selection flags are reset.

*/
/*----------------------------------------------------------------------------*/
int
uves_erase_invalid_table_rows(cpl_table *t, const char *column)
{
    int result = 0;

    assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");

    if (column == NULL)
    /* Loop through all columns */
    {
        const char *name;
        cpl_table *argument = t;
        result = 0;
        while ( (name = cpl_table_get_column_name(argument)) != NULL)
        {
            int n_deleted_rows;

            argument = NULL;
            n_deleted_rows = uves_erase_invalid_table_rows(t, name);

            if (n_deleted_rows > 0) 
            {
                uves_msg_low("%d rows with invalid '%s' removed", 
                     n_deleted_rows, name);
            }
            result += n_deleted_rows;
        }
    }
    else
    {
        assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
            "No such column: %s", column);
        
        check(( cpl_table_select_all(t),
            result = uves_table_and_selected_invalid(t, column), /* workaround here */
            cpl_table_erase_selected(t)),              /* and here */
                  "Error deleting rows");
    }
    
  cleanup:
    return result;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Extract table rows
   @param    t        Table
   @param    column   Column name
   @param    operator Logical operator
   @param    value    Value used for comparison
   @return   A new table containing the extracted rows

   A table row is extracted if and only if the value in @em column is in the
   relation @em operator to the specified @em value. The specified column
   must have type CPL_TYPE_DOUBLE or CPL_TYPE_INT or CPL_TYPE_DOUBLE.
   If integer, the integer nearest to @em value is used for the comparison.

   Also see @c cpl_table_and_selected_<type>().

*/
/*----------------------------------------------------------------------------*/
cpl_table *
uves_extract_table_rows(const cpl_table *t, const char *column,
            cpl_table_select_operator operator, double value)
{
    cpl_table *result = NULL;
    
    assure( t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure( cpl_table_has_column(t, column), CPL_ERROR_INCOMPATIBLE_INPUT,
        "No such column: %s", column);
    
    /* 1. Extract (duplicate) the entire table
       2. remove rows *not* satisfying the criterion */
    check(( result = cpl_table_duplicate(t),

        uves_select_table_rows(result, column, operator, value),
        cpl_table_not_selected(result),  /* Inverses selection */
        
        cpl_table_erase_selected(result)),
       
       "Error extracting rows");

    passure( cpl_table_count_selected(result) == cpl_table_get_nrow(result),
             "%" CPL_SIZE_FORMAT " %" CPL_SIZE_FORMAT "",
             cpl_table_count_selected(result), cpl_table_get_nrow(result) );
    
  cleanup:
    if (cpl_error_get_code() != CPL_ERROR_NONE)
    {
        uves_free_table(&result);
    }
    return result;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Sort a table by one column
   @param    t        Table
   @param    column   Column name
   @param    reverse  Flag indicating if column values are sorted descending
                      (true) or ascending (false)
   @return   CPL_ERROR_NONE iff OK

   This is a wrapper of @c cpl_table_sort().

*/
/*----------------------------------------------------------------------------*/
void
uves_sort_table_1(cpl_table *t, const char *column, bool reverse)
{
    uves_propertylist *plist = NULL;
    
    assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure(cpl_table_has_column(t, column), CPL_ERROR_ILLEGAL_INPUT, 
       "No column '%s'", column);

    check(( plist = uves_propertylist_new(),
        uves_propertylist_append_bool(plist, column, reverse)),
       "Could not create property list for sorting");

    check( uves_table_sort(t, plist), "Could not sort table");

  cleanup:
    uves_free_propertylist(&plist);
    return;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Sort a table by two columns
   @param    t        Table
   @param    column1  1st column name
   @param    column2  2nd column name
   @param    reverse1  Flag indicating if 1st column values are sorted
                       descending (true) or ascending (false)
   @param    reverse2  Flag indicating if 2nd column values are sorted
                       descending (true) or ascending (false)
   @return   CPL_ERROR_NONE iff OK

   This is a wrapper of @c cpl_table_sort(). @em column1 is the more
   significant column (i.e. values in @em column2 are compared, only if the
   values in @em column1 are equal).
*/
/*----------------------------------------------------------------------------*/
void
uves_sort_table_2(cpl_table *t, const char *column1, const char *column2, 
          bool reverse1, bool reverse2)
{
    uves_propertylist *plist = NULL;
    
    assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure(cpl_table_has_column(t, column1), CPL_ERROR_ILLEGAL_INPUT, 
       "No column '%s'", column1);
    assure(cpl_table_has_column(t, column2), CPL_ERROR_ILLEGAL_INPUT,
       "No column '%s'", column2);

    check(( plist = uves_propertylist_new(),
        uves_propertylist_append_bool(plist, column1, reverse1),
        uves_propertylist_append_bool(plist, column2, reverse2)),
       "Could not create property list for sorting");
    check( uves_table_sort(t, plist), "Could not sort table");
    
  cleanup:
    uves_free_propertylist(&plist);
    return;
}
/*----------------------------------------------------------------------------*/
/**
   @brief    Sort a table by three columns
   @param    t        Table
   @param    column1  1st column name
   @param    column2  2nd column name
   @param    column3  3rd column name
   @param    reverse1  Flag indicating if 1st column values are sorted
                       descending (true) or ascending (false)
   @param    reverse2  Flag indicating if 2nd column values are sorted
                       descending (true) or ascending (false)
   @param    reverse3  Flag indicating if 3rd column values are sorted
                       descending (true) or ascending (false)
   @return   CPL_ERROR_NONE iff OK
*/
/*----------------------------------------------------------------------------*/
void
uves_sort_table_3(cpl_table *t, const char *column1, const char *column2, 
                  const char *column3,
          bool reverse1, bool reverse2, bool reverse3)
{
    uves_propertylist *plist = NULL;
    
    assure(t != NULL, CPL_ERROR_NULL_INPUT, "Null table");
    assure(cpl_table_has_column(t, column1), CPL_ERROR_ILLEGAL_INPUT, 
       "No column '%s'", column1);
    assure(cpl_table_has_column(t, column2), CPL_ERROR_ILLEGAL_INPUT,
       "No column '%s'", column2);
    assure(cpl_table_has_column(t, column3), CPL_ERROR_ILLEGAL_INPUT,
       "No column '%s'", column3);

    check(( plist = uves_propertylist_new(),
        uves_propertylist_append_bool(plist, column1, reverse1),
        uves_propertylist_append_bool(plist, column2, reverse2),
            uves_propertylist_append_bool(plist, column3, reverse3)),
        "Could not create property list for sorting");
    check( uves_table_sort(t, plist), "Could not sort table");
    
  cleanup:
    uves_free_propertylist(&plist);
    return;
}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate memory
   @param    mem     to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free(const void *mem)
{
    cpl_free((void *)mem); /* No, it is not a bug. The cast is safe */
    return;
}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate an image and set the pointer to NULL
   @param    i        Image to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_image(cpl_image **i) 
{if(i){cpl_image_delete(*i); *i = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate an image mask and set the pointer to NULL
   @param    m        Mask to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_mask(cpl_mask **m)
{if(m){cpl_mask_delete(*m); *m = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate an image list and set the pointer to NULL
   @param    i        Image list to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_imagelist(cpl_imagelist **i)
{if(i){cpl_imagelist_delete(*i);        *i = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a table and set the pointer to NULL
   @param    t        Table to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_table(cpl_table **t)
{if(t){cpl_table_delete(*t);            *t = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a table and set the pointer to NULL
   @param    t        Table to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_table_const(const cpl_table **t)
{if(t){cpl_table_delete((cpl_table*) (*t));            *t = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a property list and set the pointer to NULL
   @param    p        Property list to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_propertylist(uves_propertylist **p)
{if(p){uves_propertylist_delete(*p);     *p = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a property list and set the pointer to NULL
   @param    p        Property list to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_propertylist_const(const uves_propertylist **p)
{if(p){uves_propertylist_delete(*p);     *p = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a property and set the pointer to NULL
   @param    p        Property to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_property(cpl_property **p)
{if(p){cpl_property_delete(*p);     *p = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a polynomial and set the pointer to NULL
   @param    p        Polynomial to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_polynomial(cpl_polynomial **p)
{if(p){cpl_polynomial_delete(*p);       *p = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a matrix and set the pointer to NULL
   @param    m        Matrix to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_matrix(cpl_matrix **m)
{if(m){cpl_matrix_delete(*m);           *m = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a parameter list and set the pointer to NULL
   @param    p        Parameter list to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_parameterlist(cpl_parameterlist **p)
{if(p){cpl_parameterlist_delete(*p);    *p = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a frameset and set the pointer to NULL
   @param    f        Frameset set to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_frameset(cpl_frameset **f)
{if(f){cpl_frameset_delete(*f);    *f = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a frame and set the pointer to NULL
   @param    f        Frame set to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_frame(cpl_frame **f)
{if(f){cpl_frame_delete(*f);    *f = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a bivector and set the pointer to NULL
   @param    b        biector to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_bivector(cpl_bivector **b)
{if(b){cpl_bivector_delete(*b);           *b = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a vector and set the pointer to NULL
   @param    v        Vector to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_vector(cpl_vector **v)
{if(v){cpl_vector_delete(*v);           *v = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a stats object and set the pointer to NULL
   @param    s        Stats object to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_stats(cpl_stats **s)
{if(s){cpl_stats_delete(*s);            *s = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Unwrap a vector and set the pointer to NULL
   @param    v        Vector to unwrap
*/
/*----------------------------------------------------------------------------*/
void uves_unwrap_vector(cpl_vector **v)
{if(v){cpl_vector_unwrap(*v);           *v = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Unwrap a vector and set the pointer to NULL
   @param    v        Vector to unwrap
*/
/*----------------------------------------------------------------------------*/
void uves_unwrap_vector_const(const cpl_vector **v)
{if(v){cpl_vector_unwrap((cpl_vector*) (*v));           *v = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Unwrap a bi-vector and set the pointer to NULL
   @param    b        Bi-vector to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_unwrap_bivector_vectors(cpl_bivector **b)
{if(b){cpl_bivector_unwrap_vectors(*b); *b = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate an array and set the pointer to NULL
   @param    a        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_array(cpl_array **a)
{if(a){cpl_array_delete(*a);           *a = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate an integer array and set the pointer to NULL
   @param    i        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_int(int **i)
{if(i){cpl_free(*i);           *i = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate an int array and set the pointer to NULL
   @param    i        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_int_const(const int **i)
{if(i){uves_free(*i);           *i = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a float array and set the pointer to NULL
   @param    f        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_float(float **f)
{if(f){cpl_free(*f);           *f = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a double array and set the pointer to NULL
   @param    i        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_double(double **d)
{if(d){cpl_free(*d);           *d = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a string array and set the pointer to NULL
   @param    i        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_string(char **s)
{if(s){cpl_free(*s);           *s = NULL;}}
/*----------------------------------------------------------------------------*/
/**
   @brief    Deallocate a double array and set the pointer to NULL
   @param    i        Array to deallocate
*/
/*----------------------------------------------------------------------------*/
void uves_free_string_const(const char **s)
{if(s){uves_free(*s);           *s = NULL;}}

/*----------------------------------------------------------------------------*/
/**
   @brief   Compute chi square
   @param   N       Number of positions
   @param   D       Dimension of x-positions
   @param   f       Function that evaluates the fit function.
   @param   a       The fit parameters.
   @param   x       Where to evaluate the fit function (N x D matrix).
   @param   y       The N values to fit.
   @param   sigma   A vector of size N containing the uncertainties of the
                    y-values. If NULL, a constant uncertainty equal to 1 is
            assumed.

   @return  chi square, or a negative number on error.

   This function calculates chi square defined as
   sum_i (y_i - f(x_i, a))^2/sigma_i^2

   Possible #_cpl_error_code_ set in this function:
   - CPL_ERROR_ILLEGAL_INPUT if the fit function could not be evaluated
*/
/*----------------------------------------------------------------------------*/

static double
get_chisq(int N, int D,
      int (*f)(const double x[], const double a[], double *result),
      const double *a,
      const double *x,
      const double *y,
      const double *sigma)
{
    double chi_sq;     /* Result */
    int i = 0;

    /* For efficiency, don't check input in this static function */

    chi_sq = 0.0;
    for (i = 0; i < N; i++)
    {
        double fx_i;
        double residual;                 /* Residual in units of uncertainty */
        const double *x_i = &(x[0+i*D]);

        /* Evaluate */
        cpl_ensure( f(x_i,
              a,
              &fx_i) == 0, CPL_ERROR_ILLEGAL_INPUT, -1.0);

        /* Accumulate */
        if (sigma == NULL)
        {
            residual = (fx_i - y[i]);
        }
        else
        {
            residual = (fx_i - y[i]) / sigma[i];
        }

        chi_sq += residual*residual;
       
    }

    return chi_sq;
}
/*----------------------------------------------------------------------------*/
/**
   @brief   Get new position in parameter space (L-M algorithm)
   @param   a       Current fit parameters.
   @param   ia      Non-NULL array defining with non-zero values which
                    parameters participate in the fit.
   @param   M       Number of fit parameters
   @param   N       Number of positions
   @param   D       Dimension of x-positions
   @param   lambda  Lambda in L-M algorithm.
   @param   f       Function that evaluates the fit function.
   @param   dfda    Function that evaluates the partial derivaties
                    of the fit function w.r.t. fit parameters.
   @param   x       The input positions (pointer to MxD matrix buffer).
   @param   y       The N values to fit.
   @param   sigma   A vector of size N containing the uncertainties of the
                    y-values. If NULL, a constant uncertainty equal to 1 is
            assumed.
   @param   partials The partial derivatives (work space).
   @param   alpha   Alpha in L-M algorithm (work space).
   @param   beta    Beta in L-M algorithm (work space).
   @param   a_da    (output) Candidate position in parameter space.

   @return  0 iff okay.

   This function computes a potentially better set of parameters @em a + @em da,
   where @em da solves the equation @em alpha(@em lambda) * @em da = @em beta .

   Possible #_cpl_error_code_ set in this function:
   - CPL_ERROR_ILLEGAL_INPUT if the fit function or its derivative could
   not be evaluated.
   - CPL_ERROR_SINGULAR_MATRIX if @em alpha is singular.

*/
/*----------------------------------------------------------------------------*/
static int
get_candidate(const double *a, const int ia[],
          int M, int N, int D,
          double lambda,
          int    (*f)(const double x[], const double a[], double *result),
          int (*dfda)(const double x[], const double a[], double result[]),
          const double *x,
          const double *y,
          const double *sigma,
          double *partials,
          cpl_matrix *alpha,
          cpl_matrix *beta,
          double *a_da)
{
    int Mfit = 0;         /* Number of non-constant fit parameters */
    cpl_matrix *da;       /* Solution of   alpha * da = beta */
    double *alpha_data;
    double *beta_data;
    double *da_data;
    int i, imfit = 0;
    int j, jmfit = 0;
    int k = 0;

    /* For efficiency, don't check input in this static function */

    Mfit = cpl_matrix_get_nrow(alpha);

    alpha_data    = cpl_matrix_get_data(alpha);
    beta_data     = cpl_matrix_get_data(beta);
   
    /* Build alpha, beta:
     *
     *  alpha[i,j] = sum_{k=1,N} (sigma_k)^-2 * df/da_i * df/da_j  *
     *                           (1 + delta_ij lambda) ,
     *
     *   beta[i]   = sum_{k=1,N} (sigma_k)^-2 * ( y_k - f(x_k) ) * df/da_i
     *
     * where (i,j) loop over the non-constant parameters (0 to Mfit-1),
     * delta is Kronecker's delta, and all df/da are evaluated in x_k
     */

    cpl_matrix_fill(alpha, 0.0);
    cpl_matrix_fill(beta , 0.0);

    for (k = 0; k < N; k++)
    {
        double sm2 = 0.0;                /* (sigma_k)^-2 */
        double fx_k = 0.0;               /* f(x_k)       */
        const double *x_k = &(x[0+k*D]); /* x_k          */

        if (sigma == NULL)
        {
            sm2 = 1.0;
        }
        else
        {
            sm2 = 1.0 / (sigma[k] * sigma[k]);
        }
        
        /* Evaluate f(x_k) */
        cpl_ensure( f(x_k, a, &fx_k) == 0, CPL_ERROR_ILLEGAL_INPUT, -1);

        /* Evaluate (all) df/da (x_k) */
        cpl_ensure( dfda(x_k, a, partials) == 0, 
            CPL_ERROR_ILLEGAL_INPUT, -1);

        for (i = 0, imfit = 0; i < M; i++)
        {
            if (ia[i] != 0)
            {
                /* Beta */
                beta_data[imfit] +=
                sm2 * (y[k] - fx_k) * partials[i];
                
                /* Alpha is symmetrical, so compute
                   only lower-left part */
                for (j = 0, jmfit = 0; j < i; j++)
                {
                    if (ia[j] != 0)
                    {
                        alpha_data[jmfit + imfit*Mfit] +=
                        sm2 * partials[i] * 
                        partials[j];
                        
                        jmfit += 1;
                    }
                }
                
                /* Alpha, diagonal terms */
                j = i;
                jmfit = imfit;
                
                alpha_data[jmfit + imfit*Mfit] += 
                sm2 * partials[i] *
                partials[j] * (1 + lambda);
                
                imfit += 1;
            }
        }

        assert( imfit == Mfit );
    }

    /* Create upper-right part of alpha */
    for (i = 0, imfit = 0; i < M; i++)
    {
        if (ia[i] != 0)
        {
            for (j = i+1, jmfit = imfit+1; j < M; j++)
            {
                if (ia[j] != 0)
                {
                    alpha_data[jmfit + imfit*Mfit] = 
                    alpha_data[imfit + jmfit*Mfit];
                    
                    jmfit += 1;
                }
            }
            assert( jmfit == Mfit );

            imfit += 1;
        }
    }
    assert( imfit == Mfit );
    
    da = cpl_matrix_solve(alpha, beta);

    cpl_ensure(da != NULL, cpl_error_get_code(), -1);

    /* Create a+da vector by adding a and da */
    da_data   = cpl_matrix_get_data(da);

    for (i = 0, imfit = 0; i < M; i++)
    {
        if (ia[i] != 0)
        {
            a_da[i] = a[i] + da_data[0 + imfit*1];
            
            imfit += 1;
        }
        else
        {
            a_da[i] = a[i];
        }
    }
    
    assert( imfit == Mfit );

    cpl_matrix_delete(da);

    return 0;
}

/**@}*/