GSL 0.4.23

A rust binding for the GSL (the GNU scientific library)
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
//
// A rust binding for the GSL library by Guillaume Gomez (guillaume1.gomez@gmail.com)
//

/*!
#Nonlinear Least-Squares Fitting

This chapter describes functions for multidimensional nonlinear least-squares fitting. The library provides low level components for a 
variety of iterative solvers and convergence tests. These can be combined by the user to achieve the desired solution, with full access to 
the intermediate steps of the iteration. Each class of methods uses the same framework, so that you can switch between solvers at runtime 
without needing to recompile your program. Each instance of a solver keeps track of its own state, allowing the solvers to be used in 
multi-threaded programs.

##Overview

The problem of multidimensional nonlinear least-squares fitting requires the minimization of the squared residuals of n functions, f_i, in 
p parameters, x_i,

\Phi(x) = (1/2) || F(x) ||^2
        = (1/2) \sum_{i=1}^{n} f_i(x_1, ..., x_p)^2 
All algorithms proceed from an initial guess using the linearization,

\psi(p) = || F(x+p) || ~=~ || F(x) + J p ||
where x is the initial point, p is the proposed step and J is the Jacobian matrix J_{ij} = d f_i / d x_j. Additional strategies are used to 
enlarge the region of convergence. These include requiring a decrease in the norm ||F|| on each step or using a trust region to avoid steps 
which fall outside the linear regime.

To perform a weighted least-squares fit of a nonlinear model Y(x,t) to data (t_i, y_i) with independent Gaussian errors \sigma_i, use 
function components of the following form,

f_i = (Y(x, t_i) - y_i) / \sigma_i
Note that the model parameters are denoted by x in this chapter since the non-linear least-squares algorithms are described geometrically 
(i.e. finding the minimum of a surface). The independent variable of any data to be fitted is denoted by t.

With the definition above the Jacobian is J_{ij} =(1 / \sigma_i) d Y_i / d x_j, where Y_i = Y(x,t_i).

##High Level Driver

These routines provide a high level wrapper that combine the iteration and convergence testing for easy use.
!*/

use ffi;
use num::Float;
use libc::c_void;

pub struct MultiFitFunction<'r, T:'r> {
    pub f: fn(x: &::VectorF64, params: &mut T, f: &::VectorF64) -> ::Value,
    /// number of functions
    pub n: u64,
    /// number of independent variables
    pub p: u64,
    pub params: &'r mut T
}

pub struct MultiFitFdfSolver<'r, T:'r> {
    _type: &'r MultiFitFdfSolverType<T>,
    fdf: *mut c_void, //MultiFitFunctionFdf<'r, T>,
    pub x: ::VectorF64,
    pub f: ::VectorF64,
    pub j: ::MatrixF64,
    pub dx: ::VectorF64,
    state: LmderStateT
}

impl<'r, T> MultiFitFdfSolver<'r, T> {
    /// This function returns a pointer to a newly allocated instance of a solver of type T for n observations and p parameters. The number
    /// of observations n must be greater than or equal to parameters p.
    pub fn new(_type: &'r MultiFitFdfSolverType<T>, n: u64, p: u64) -> Option<MultiFitFdfSolver<'r, T>> {
        let mut r = MultiFitFdfSolver {
            _type: _type,
            fdf: ::std::ptr::null_mut(),
            x: ::VectorF64::new(p).unwrap(),
            f: ::VectorF64::new(n).unwrap(),
            j: ::MatrixF64::new(n, p).unwrap(),
            dx: ::VectorF64::new(p).unwrap(),
            state: LmderStateT::new(n, p)
        };
        if ((*_type).alloc)(&mut r.state, n, p) == ::Value::Success {
            Some(r)
        } else {
            None
        }
    }

    /// This function initializes, or reinitializes, an existing solver s to use the function f and the initial guess x.
    pub fn set(&mut self, f: &mut MultiFitFunctionFdf<'r, T>, x: &::VectorF64) -> ::Value {
        if self.f.len() != f.n {
            rgsl_error!("function size does not match solver", ::Value::BadLen);
            return ::Value::BadLen;
        }

        if self.x.len() != x.len() {
            rgsl_error!("vector length does not match solver", ::Value::BadLen);
            return ::Value::BadLen;
        }  

        self.fdf = unsafe { ::std::mem::transmute(f) };
        self.x.copy_from(x);

        unsafe {
            (self._type.set)(&mut self.state, ::std::mem::transmute(self.fdf), &mut self.x, &mut self.f, &mut self.j, &mut self.dx)
        }
    }

    pub fn name(&self) -> &'static str {
        self._type.name
    }

    /// This function performs a single iteration of the solver s. If the iteration encounters an unexpected problem then an error code
    /// will be returned. The solver maintains a current estimate of the best-fit parameters at all times.
    pub fn iterate(&mut self) -> ::Value {
        unsafe {
            (self._type.iterate)(&mut self.state, ::std::mem::transmute(self.fdf), &mut self.x, &mut self.f, &mut self.j, &mut self.dx)
        }
    }

    /// This function returns the current position (i.e. best-fit parameters) s->x of the solver s.
    pub fn position(&'r self) -> &'r ::VectorF64 {
        &self.x
    }

    /// These functions iterate the solver s for a maximum of maxiter iterations. After each iteration, the system is tested for convergence
    /// using gsl_multifit_test_delta with the error tolerances epsabs and epsrel.
    #[allow(unused_assignments)]
    pub fn driver(&mut self, max_iter: u64, epsabs: f64, epsrel: f64) -> ::Value {
        let mut status = ::Value::Success;
        let mut iter = 0u64;

        loop {
            status = self.iterate();

            if status != ::Value::Success {
                break;
            }

            /* test for convergence */
            status = gsl_multifit_test_delta(&self.dx, &self.x, epsabs, epsrel);
            iter += 1;
            if status != ::Value::Continue || iter >= max_iter {
                break;
            }
        }

        status
    }
}

impl<'r, T> Drop for MultiFitFdfSolver<'r, T> {
    fn drop(&mut self) {
        //self.s = ::std::ptr::null_mut();
    }
}

#[allow(dead_code)]
pub struct MultiFitFdfSolverType<T> {
    name: &'static str,
    //size: u64,
    alloc: fn(state: &mut LmderStateT, n: u64, p: u64) -> ::Value,
    set: fn(state: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, j: &mut ::MatrixF64,
        dx: &mut ::VectorF64) -> ::Value,
    iterate: fn(state: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64,
        j: &mut ::MatrixF64, dx: &mut ::VectorF64) -> ::Value,
    free: fn(state: &mut LmderStateT)
}

impl<T> MultiFitFdfSolverType<T> {
    pub fn lmder() -> MultiFitFdfSolverType<T> {
        MultiFitFdfSolverType {
            name: "lmder",
            alloc: lmder_alloc,
            set: lmder_set,
            iterate: lmder_iterate,
            free: lmder_free
        }
    }

    pub fn lmsder() -> MultiFitFdfSolverType<T> {
        MultiFitFdfSolverType {
            name: "lmsder",
            alloc: lmder_alloc,
            set: lmsder_set,
            iterate: lmder_iterate,
            free: lmder_free
        }
    }
}

pub struct MultiFitFunctionFdf<'r, T:'r> {
    pub f: fn(x: &::VectorF64, params: &mut T, f: &mut ::VectorF64) -> ::Value,
    pub df: Option<fn(x: &::VectorF64, params: &mut T, df: &mut ::MatrixF64) -> ::Value>,
    pub fdf: Option<fn(x: &::VectorF64, params: &mut T, f: &mut ::VectorF64, df: &mut ::MatrixF64) -> ::Value>,
    pub n: u64,
    pub p: u64,
    pub params: &'r mut T
}

struct LmderStateT {
    iter: u64,
    xnorm: f64,
    fnorm: f64,
    delta: f64,
    par: f64,
    r: ::MatrixF64,
    tau: ::VectorF64,
    diag: ::VectorF64,
    qtf: ::VectorF64,
    newton: ::VectorF64,
    gradient: ::VectorF64,
    x_trial: ::VectorF64,
    f_trial: ::VectorF64,
    df: ::VectorF64,
    sdiag: ::VectorF64,
    rptdx: ::VectorF64,
    w: ::VectorF64,
    work1: ::VectorF64,
    perm: ::Permutation
}

impl LmderStateT {
    fn new(n: u64, p: u64) -> LmderStateT {
        LmderStateT {
            iter: 0u64,
            xnorm: 0f64,
            fnorm: 0f64,
            delta: 0f64,
            par: 0f64,
            r: ::MatrixF64::new(n, p).unwrap(),
            tau: ::VectorF64::new(if p < n {p} else {n}).unwrap(),
            diag: ::VectorF64::new(p).unwrap(),
            qtf: ::VectorF64::new(n).unwrap(),
            newton: ::VectorF64::new(p).unwrap(),
            gradient: ::VectorF64::new(p).unwrap(),
            x_trial: ::VectorF64::new(p).unwrap(),
            f_trial: ::VectorF64::new(n).unwrap(),
            df: ::VectorF64::new(n).unwrap(),
            sdiag: ::VectorF64::new(p).unwrap(),
            rptdx: ::VectorF64::new(n).unwrap(),
            w: ::VectorF64::new(n).unwrap(),
            work1: ::VectorF64::new(p).unwrap(),
            perm: ::Permutation::new(p).unwrap()
        }
    }
}

#[allow(unused_variables)]
// I'm not sure if I'll keep this function or not...
fn lmder_alloc(state: &mut LmderStateT, n: u64, p: u64) -> ::Value {
    /*state.r = match ::MatrixF64::new(n, p) {
        Some(m) => m,
        None => {
            rgsl_error!("failed to allocate space for r", ::Value::NoMem);
            ::MatrixF64::new(0, 0).unwrap() // to allow compilation
        }
    };

    state.tau = match ::VectorF64::new(if n < p {n} else {p}) {
        Some(t) => t,
        None => {
            rgsl_error!("failed to allocate space for tau", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.diag = match ::VectorF64::new(p) {
        Some(d) => d,
        None => {
            rgsl_error!("failed to allocate space for diag", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.qtf = match ::VectorF64::new(n) {
        Some(q) => q,
        None => {
            rgsl_error!("failed to allocate space for qtf", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.newton = match ::VectorF64::new(p) {
        Some(n) => n,
        None => {
            rgsl_error!("failed to allocate space for newton", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.gradient = match ::VectorF64::new(p) {
        Some(g) => g,
        None => {
            rgsl_error!("failed to allocate space for gradient", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.x_trial = match ::VectorF64::new(p) {
        Some(x) => x,
        None => {
            rgsl_error!("failed to allocate space for x_trial", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.f_trial = match ::VectorF64::new(n) {
        Some(f) => f,
        None => {
            rgsl_error!("failed to allocate space for f_trial", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.df = match ::VectorF64::new(n) {
        Some(d) => d,
        None => {
            rgsl_error!("failed to allocate space for df", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.sdiag = match ::VectorF64::new(p) {
        Some(s) => s,
        None => {
            rgsl_error!("failed to allocate space for sdiag", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.rptdx = match ::VectorF64::new(n) {
        Some(s) => s,
        None => {
            rgsl_error!("failed to allocate space for rptdx", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.w = match ::VectorF64::new(n) {
        Some(w) => w,
        None => {
            rgsl_error!("failed to allocate space for w", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.work1 = match ::VectorF64::new(p) {
        Some(w) => w,
        None => {
            rgsl_error!("failed to allocate space for work1", ::Value::NoMem);
            ::VectorF64::new(0).unwrap() // to allow compilation
        }
    };

    state.perm = match ::Permutation::new(p) {
        Some(p) => p,
        None => {
            rgsl_error!("failed to allocate space for perm", ::Value::NoMem);
            ::Permutation::new(0).unwrap() // to allow compilation
        }
    };*/

    ::Value::Success
}

#[allow(unused_variables)]
fn lmder_free(state: &mut LmderStateT) {
    /*::std::mem::drop(state.perm);
    ::std::mem::drop(state.work1);
    ::std::mem::drop(state.w);
    ::std::mem::drop(state.rptdx);
    ::std::mem::drop(state.sdiag);
    ::std::mem::drop(state.df);
    ::std::mem::drop(state.f_trial);
    ::std::mem::drop(state.x_trial);
    ::std::mem::drop(state.gradient);
    ::std::mem::drop(state.newton);
    ::std::mem::drop(state.qtf);
    ::std::mem::drop(state.diag);
    ::std::mem::drop(state.tau);
    ::std::mem::drop(state.r);*/
}

fn lmder_set<T>(vstate: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, J: &mut ::MatrixF64,
    dx: &mut ::VectorF64) -> ::Value {
    set(vstate, fdf, x, f, J, dx, 0)
}

fn lmsder_set<T>(vstate: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, J: &mut ::MatrixF64,
    dx: &mut ::VectorF64) -> ::Value {
    set(vstate, fdf, x, f, J, dx, 1)
}

fn lmder_iterate<T>(vstate: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, J: &mut ::MatrixF64,
    dx: &mut ::VectorF64) -> ::Value {
    iterate(vstate, fdf, x, f, J, dx, 0)
}

#[allow(dead_code)]
fn lmsder_iterate<T>(vstate: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, J: &mut ::MatrixF64,
    dx: &mut ::VectorF64) -> ::Value {
    iterate(vstate, fdf, x, f, J, dx, 1)
}

fn compute_diag(J: &::MatrixF64, diag: &mut ::VectorF64) {
    let n = J.size1();
    let p = J.size2();

    for j in 0..p {
        let mut sum = 0f64;

        for i in 0..n {
            let jij = J.get(i, j);

            sum += jij * jij;
        }
        if sum == 0f64 {
            sum = 0f64;
        }

        diag.set(j, sum.sqrt());
    }
}

fn update_diag(J: &::MatrixF64, diag: &mut ::VectorF64) {
    let n = diag.len();

    for j in 0..n {
        let mut sum = 0f64;

        for i in 0..n {
          let jij = J.get(i, j);
          
          sum += jij * jij;
        }
        if sum == 0f64 {
            sum = 1f64;
        }

        let cnorm = sum.sqrt();
        let diagj = diag.get(j);

        if cnorm > diagj {
            diag.set(j, cnorm);
        }
    }
}

fn scaled_enorm(d: &::VectorF64, f: &::VectorF64) -> f64 {
    let mut e2 = 0f64;
    let n = f.len();

    for i in 0..n {
        let fi = f.get(i);
        let di = d.get(i);
        let u = di * fi;

        e2 += u * u;
    }
    e2.sqrt()
}

fn compute_delta(diag: &mut ::VectorF64, x: &mut ::VectorF64) -> f64 {
    let Dx = scaled_enorm(diag, x);
    let factor = 100f64;  /* generally recommended value from MINPACK */

    if Dx > 0f64 {
        factor * Dx
    } else {
        factor
    }
}

fn set<T>(state: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, j: &mut ::MatrixF64,
    dx: &mut ::VectorF64, scale: i32) -> ::Value {
    let mut signum = 0i32;

    /* Evaluate function at x */
    /* return immediately if evaluation raised error */
    {
        let status = match fdf.fdf {
            Some(func) => {
                func(x, fdf.params, f, j)
            }
            None => {
                /* finite difference approximation */
                gsl_multifit_fdfsolver_dif_fdf(x, fdf, f, j)
            }
        };

        if status != ::Value::Success {
            return status;
        }
    }

    state.par = 0f64;
    state.iter = 1;
    state.fnorm = enorm(f);

    dx.set_all(0f64);

    /* store column norms in diag */
    if scale != 0 {
        compute_diag(j, &mut state.diag);
    } else {
        state.diag.set_all(1f64);
    }

    /* set delta to 100 |D x| or to 100 if |D x| is zero */
    state.xnorm = scaled_enorm(&state.diag, x);
    state.delta = compute_delta(&mut state.diag, x);

    /* Factorize J into QR decomposition */
    state.r.copy_from(j);
    ::linear_algebra::QRPT_decomp(&state.r, &state.tau, &state.perm, &mut signum, &state.work1);

    state.rptdx.set_zero();
    state.w.set_zero();

    /* Zero the trial vector, as in the alloc function */

    state.f_trial.set_zero();

    /*#ifdef DEBUG
    printf("r = "); gsl_matrix_fprintf(stdout, r, "%g");
    printf("perm = "); gsl_permutation_fprintf(stdout, perm, "%d");
    printf("tau = "); gsl_vector_fprintf(stdout, tau, "%g");
    #endif*/

    ::Value::Success
}

fn enorm(f: &::VectorF64) -> f64 {
    ::blas::level1::dnrm2(f)
}

fn compute_gradient_direction(r: &::MatrixF64, p: &::Permutation, qtf: &::VectorF64, diag: &::VectorF64, g: &mut ::VectorF64) {
    let n = r.size2();

    for j in 0..n {
        let mut sum = 0f64;

        for i in 0..(j + 1) {
            sum += r.get(i, j) * qtf.get(i);
        }

        {
            let pj = p.get(j);
            let dpj = diag.get(pj);

            g.set(j, sum / dpj);
        }
    }
}

fn compute_trial_step(x: &mut ::VectorF64, dx: &mut ::VectorF64, x_trial: &mut ::VectorF64) {
    let n = x.len();

    for i in 0..n {
        let pi = dx.get(i);
        let xi = x.get(i);

        x_trial.set(i, xi + pi);
    }
}

fn compute_actual_reduction(fnorm: f64, fnorm1: f64) -> f64 {
    if 0.1f64 * fnorm1 < fnorm {
        let u = fnorm1 / fnorm;

        1f64 - u * u
    } else {
        -1f64
    }
}

fn compute_rptdx(r: &::MatrixF64, p: &::Permutation, dx: &::VectorF64, rptdx: &mut ::VectorF64) {
    let n = dx.len();

    for i in 0..n {
        let mut sum = 0f64;

        for j in i..n {
            let pj = p.get(j);

            sum += r.get(i, j) * dx.get(pj);
        }

        rptdx.set(i, sum);
    }
}

#[allow(unused_assignments)]
fn iterate<T>(state: &mut LmderStateT, fdf: &mut MultiFitFunctionFdf<T>, x: &mut ::VectorF64, f: &mut ::VectorF64, J: &mut ::MatrixF64,
    dx: &mut ::VectorF64, scale: i32) -> ::Value {
    let mut prered = 0f64;
    let mut actred = 0f64;
    let mut pnorm = 0f64;
    let mut fnorm1 = 0f64;
    let mut fnorm1p = 0f64;
    let mut gnorm = 0f64;
    let mut dirder = 0f64;

    let mut iter = 0i32;

    let p1 = 0.1f64;
    let p25 = 0.25f64;
    let p5 = 0.5f64;
    let p75 = 0.75f64;
    let p0001 = 0.0001f64;

    if state.fnorm == 0f64 {
        return ::Value::Success;
    }

    /* Compute qtf = Q^T f */

    state.qtf.copy_from(f);
    ::linear_algebra::QR_QTvec(&state.r, &state.tau, &state.qtf);

    /* Compute norm of scaled gradient */

    compute_gradient_direction(&state.r, &state.perm, &state.qtf, &state.diag, &mut state.gradient);

    { 
        let iamax = ::blas::level1::idamax(&state.gradient);

        gnorm = unsafe { (state.gradient.get(iamax as u64) / state.fnorm).abs() };
    }

    /* Determine the Levenberg-Marquardt parameter */

    loop {
        iter += 1;

        {
            let status = lmpar(&mut state.r, &state.perm, &state.qtf, &state.diag, state.delta, &mut (state.par), &mut state.newton,
                &mut state.gradient, &mut state.sdiag, dx, &mut state.w);

            if status != ::Value::Success {
                return status;
            }
        }

        /* Take a trial step */

        dx.scale(-1f64); /* reverse the step to go downhill */

        compute_trial_step(x, dx, &mut state.x_trial);

        pnorm = scaled_enorm(&state.diag, dx);

        if state.iter == 1 {
            if pnorm < state.delta {
        /*#ifdef DEBUG
              printf("set delta = pnorm = %g\n" , pnorm);
        #endif*/
                state.delta = pnorm;
            }
        }

        /* Evaluate function at x + p */
        /* return immediately if evaluation raised error */
        {
            let status = ((*fdf).f)(&state.x_trial, fdf.params, &mut state.f_trial); //GSL_MULTIFIT_FN_EVAL_F (fdf, x_trial, f_trial);
            
            if status != ::Value::Success {
                return status;
            }
        }

        fnorm1 = enorm(&state.f_trial);

        /* Compute the scaled actual reduction */

        actred = compute_actual_reduction(state.fnorm, fnorm1);

        /*#ifdef DEBUG
            printf("lmiterate: fnorm = %g fnorm1 = %g  actred = %g\n", state->fnorm, fnorm1, actred);
            printf("r = "); gsl_matrix_fprintf(stdout, r, "%g");
            printf("perm = "); gsl_permutation_fprintf(stdout, perm, "%d");
            printf("dx = "); gsl_vector_fprintf(stdout, dx, "%g");
        #endif*/

        /* Compute rptdx = R P^T dx, noting that |J dx| = |R P^T dx| */

        compute_rptdx(&state.r, &state.perm, dx, &mut state.rptdx);

        /*#ifdef DEBUG
        printf("rptdx = "); gsl_vector_fprintf(stdout, rptdx, "%g");
        #endif*/

        fnorm1p = enorm(&state.rptdx);

        /* Compute the scaled predicted reduction = |J dx|^2 + 2 par |D dx|^2 */

        { 
            let t1 = fnorm1p / state.fnorm;
            let t2 = (state.par.sqrt() * pnorm) / state.fnorm;

            prered = t1 * t1 + t2 * t2 / p5;
            dirder = -(t1 * t1 + t2 * t2);
        }

        /* compute the ratio of the actual to predicted reduction */

        let ratio = if prered > 0f64 {
            actred / prered
        } else {
            0f64
        };

        /*#ifdef DEBUG
        printf("lmiterate: prered = %g dirder = %g ratio = %g\n", prered, dirder,ratio);
        #endif*/

        /* update the step bound */
        if ratio > p25 {
    /*#ifdef DEBUG
          printf("ratio > p25\n");
    #endif*/
            if state.par == 0f64 || ratio >= p75 {
                state.delta = pnorm / p5;
                state.par *= p5;
    /*#ifdef DEBUG
              printf("updated step bounds: delta = %g, par = %g\n", state->delta, state->par);
    #endif*/
            }
        } else {
            let mut temp = if actred >= 0f64 {p5} else {p5 * dirder / (dirder + p5 * actred)};

    /*#ifdef DEBUG
          printf("ratio < p25\n");
    #endif*/

            if p1 * fnorm1 >= state.fnorm || temp < p1 {
                temp = p1;
            }

            state.delta = temp * state.delta.min(pnorm / p1);

            state.par /= temp;
    /*#ifdef DEBUG
          printf("updated step bounds: delta = %g, par = %g\n", state->delta, state->par);
    #endif*/
        }


        /* test for successful iteration, termination and stringent tolerances */

        if ratio >= p0001 {
            x.copy_from(&state.x_trial);
            f.copy_from(&state.f_trial);

            /* return immediately if evaluation raised error */
            {
                let status = match fdf.df {
                    Some(func) => {
                        func(&state.x_trial, fdf.params, J)
                    }
                    None => {
                        gsl_multifit_fdfsolver_dif_df(&mut state.x_trial, fdf, &mut state.f_trial, J)
                    }
                };

                if status != ::Value::Success {
                    return status;
                }
            }

            /* wa2_j  = diag_j * x_j */
            state.xnorm = scaled_enorm(&state.diag, x);
            state.fnorm = fnorm1;
            state.iter += 1;

            /* Rescale if necessary */
            if scale != 0 {
                update_diag(J, &mut state.diag);
            }

            {
                let mut signum = 0;

                state.r.copy_from(J);
                ::linear_algebra::QRPT_decomp(&state.r, &state.tau, &state.perm, &mut signum, &state.work1);
            }

            return ::Value::Success;
        } else if actred.abs() <= ::DBL_EPSILON  && prered <= ::DBL_EPSILON  && p5 * ratio <= 1f64 {
            return ::Value::TolF;
        } else if state.delta <= ::DBL_EPSILON * state.xnorm {
            return ::Value::TolX;
        } else if gnorm <= ::DBL_EPSILON {
            return ::Value::TolG;
        } else if iter >= 10 {
            /* stop inner loop if successful */
            break;
        }
    }

    return ::Value::NoProg;
}

fn gsl_multifit_test_delta(dx: &::VectorF64, x: &::VectorF64, epsabs: f64, epsrel: f64) -> ::Value {
    let mut ok = false;

    if epsrel < 0f64 {
        rgsl_error!("relative tolerance is negative", ::Value::BadTol);
        return ::Value::BadTol;
    }

    for i in 0..x.len() {
        let xi = x.get(i);
        let dxi = dx.get(i);
        let tolerance = epsabs + epsrel * xi.abs();

        if dxi.abs() < tolerance {
            ok = true;
        } else {
            ok = false;
            break;
        }
    }

    if ok {
        ::Value::Success
    } else {
        ::Value::Continue
    }
}

/*
gsl_multifit_fdfsolver_dif_fdf()
  Compute function values (analytic) and approximate Jacobian using finite
differences
Inputs: x      - parameter vector
        fdf    - fdf
        f      - (output) function values f_i(x)
        J      - (output) approximate Jacobian matrix
Return: success or error
*/

fn gsl_multifit_fdfsolver_dif_fdf<T>(x: &mut ::VectorF64, fdf: &mut MultiFitFunctionFdf<T>, f: &mut ::VectorF64,
    j: &mut ::MatrixF64) -> ::Value {
    let mut status = ((*fdf).f)(x, fdf.params, f); // GSL_MULTIFIT_FN_EVAL_F(fdf, x, f);
    if status != ::Value::Success {
        return status;
    }

    status = fdjac(x, fdf, f, j);
    if status != ::Value::Success {
        return status;
    }

    status
}

/*
fdjac()
  Compute approximate Jacobian using forward differences
Inputs: x   - parameter vector
        fdf - fdf struct
        f   - (input) vector of function values f_i(x)
        J   - (output) Jacobian matrix
Return: success or error
*/

fn fdjac<T>(x: &mut ::VectorF64, fdf: &mut MultiFitFunctionFdf<T>, f: &::VectorF64, jm: &mut ::MatrixF64) -> ::Value {
    let mut status = ::Value::Success;
    let epsfcn = 0f64;
    let eps = (epsfcn.max(::DBL_EPSILON)).sqrt();

    for j in 0..fdf.p {
        let xj = x.get(j);

        /* use column j of J as temporary storage for f(x + dx) */
        let (mut v, _) = jm.get_col(j).unwrap();

        let mut h = eps * xj.abs();
        if h == 0f64 {
            h = eps;
        }

        /* perturb x_j to compute forward difference */
        x.set(j, xj + h);

        status = ((*fdf).f)(x, fdf.params, &mut v); //GSL_MULTIFIT_FN_EVAL_F (fdf, x, &v.vector);
        if status != ::Value::Success {
            return status;
        }

        /* restore x_j */
        x.set(j, xj);

        h = 1f64 / h;
        for i in 0..fdf.n {
            let fnext = v.get(i);
            let fi = f.get(i);

            jm.set(i, j, (fnext - fi) * h);
        }
    }

    status
}

/*
gsl_multifit_fdfsolver_dif_df()
  Compute approximate Jacobian using finite differences
Inputs: x   - parameter vector
        fdf - fdf
        f   - (input) function values f_i(x)
        J   - (output) approximate Jacobian matrix
Return: success or error
*/

fn gsl_multifit_fdfsolver_dif_df<T>(x: &mut ::VectorF64, fdf: &mut MultiFitFunctionFdf<T>, f: &::VectorF64,
    J: &mut ::MatrixF64) -> ::Value {
    fdjac(x, fdf, f, J)
}

#[allow(unused_assignments)]
fn lmpar(r: &mut ::MatrixF64, perm: &::Permutation, qtf: &::VectorF64, diag: &::VectorF64, delta: f64, par_inout: &mut f64,
    newton: &mut ::VectorF64, gradient: &mut ::VectorF64, sdiag: &mut ::VectorF64, x: &mut ::VectorF64,
    w: &mut ::VectorF64) -> ::Value {
    let mut dxnorm = 0f64;
    let mut gnorm = 0f64;
    let mut fp = 0f64;
    let mut fp_old = 0f64;
    let mut par_lower = 0f64;
    let mut par_upper = 0f64;
    let mut par_c = 0f64;

    let mut par = *par_inout;
    let mut iter = 0u64;

    /*#ifdef DEBUG
    printf("ENTERING lmpar\n");
    #endif*/


    compute_newton_direction(r, perm, qtf, newton);

    /*#ifdef DEBUG
    printf ("newton = ");
    gsl_vector_fprintf (stdout, newton, "%g");
    printf ("\n");

    printf ("diag = ");
    gsl_vector_fprintf (stdout, diag, "%g");
    printf ("\n");
    #endif*/

    /* Evaluate the function at the origin and test for acceptance of
     the Gauss-Newton direction. */

    dxnorm = scaled_enorm(diag, newton);

    fp = dxnorm - delta;

    /*#ifdef DEBUG
    printf ("dxnorm = %g, delta = %g, fp = %g\n", dxnorm, delta, fp);
    #endif*/

    if fp <= 0.1 * delta {
        x.copy_from(newton);
        /*#ifdef DEBUG
          printf ("took newton (fp = %g, delta = %g)\n", fp, delta);
        #endif*/

        *par_inout = 0f64;

        return ::Value::Success;
    }

    /*#ifdef DEBUG
    printf ("r = ");
    gsl_matrix_fprintf (stdout, r, "%g");
    printf ("\n");

    printf ("newton = ");
    gsl_vector_fprintf (stdout, newton, "%g");
    printf ("\n");

    printf ("dxnorm = %g\n", dxnorm);
    #endif*/

    compute_newton_bound(r, newton, dxnorm, perm, diag, w);

    /*#ifdef DEBUG
    printf("perm = "); gsl_permutation_fprintf(stdout, perm, "%d");

    printf ("diag = ");
    gsl_vector_fprintf (stdout, diag, "%g");
    printf ("\n");

    printf ("w = ");
    gsl_vector_fprintf (stdout, w, "%g");
    printf ("\n");
    #endif*/


    {
        let wnorm = enorm(w);
        let phider = wnorm * wnorm;

        /* w == zero if r rank-deficient, 
           then set lower bound to zero form MINPACK, lmder.f 
           Hans E. Plesser 2002-02-25 (hans.plesser@itf.nlh.no) */
        par_lower = if wnorm > 0f64 {
            fp / (delta * phider)
        } else {
            0f64
        };
    }

    /*#ifdef DEBUG
    printf("par       = %g\n", par      );
    printf("par_lower = %g\n", par_lower);
    #endif*/

    compute_gradient_direction(r, perm, qtf, diag, gradient);

    gnorm = enorm(gradient);

    /*#ifdef DEBUG
    printf("gradient = "); gsl_vector_fprintf(stdout, gradient, "%g"); printf("\n");
    printf("gnorm = %g\n", gnorm);
    #endif*/

    par_upper =  gnorm / delta;

    if par_upper == 0f64 {
        par_upper = ::DBL_MIN / delta.min(0.1f64);
    }

    /*#ifdef DEBUG
    printf("par_upper = %g\n", par_upper);
    #endif*/

    if par > par_upper {
        /*#ifdef DEBUG
        printf("set par to par_upper\n");
        #endif*/

        par = par_upper;
    } else if par < par_lower {
        /*#ifdef DEBUG
        printf("set par to par_lower\n");
        #endif*/

        par = par_lower;
    }

    if par == 0f64 {
        par = gnorm / dxnorm;
        /*#ifdef DEBUG
        printf("set par to gnorm/dxnorm = %g\n", par);
        #endif*/
    }

    /* Beginning of iteration */
    loop {
        iter += 1;

        /*#ifdef DEBUG
        printf("lmpar iteration = %d\n", iter);
        #endif*/

        //#ifdef BRIANSFIX
        /* Seems like this is described in the paper but not in the MINPACK code */

        // I keep the BRIANSFIX for the moment
        if par < par_lower || par > par_upper {
            par = (par_lower * par_upper).sqrt().max(0.001f64 * par_upper);
        }
        //#endif

        /* Evaluate the function at the current value of par */

        if par == 0f64 {
            par = (0.001f64 * par_upper).max(::DBL_MIN);
            /*#ifdef DEBUG
              printf("par = 0, set par to  = %g\n", par);
            #endif*/
        }

        /* Compute the least squares solution of [ R P x - Q^T f, sqrt(par) D x]
         for A = Q R P^T */

        /*#ifdef DEBUG
        printf ("calling qrsolv with par = %g\n", par);
        #endif*/

        {
            let sqrt_par = par.sqrt();

            qrsolv(r, perm, sqrt_par, diag, qtf, x, sdiag, w);
        }

        dxnorm = scaled_enorm(diag, x);

        fp_old = fp;

        fp = dxnorm - delta;

        /*#ifdef DEBUG
        printf ("After qrsolv dxnorm = %g, delta = %g, fp = %g\n", dxnorm, delta, fp);
        printf ("sdiag = ") ; gsl_vector_fprintf(stdout, sdiag, "%g"); printf("\n");
        printf ("x = ") ; gsl_vector_fprintf(stdout, x, "%g"); printf("\n");
        printf ("r = ") ; gsl_matrix_fprintf(stdout, r, "%g"); printf("\nXXX\n");
        #endif*/

        /* If the function is small enough, accept the current value of par */

        if fp.abs() <= 0.1f64 * delta {
            break;
        }

        if par_lower == 0f64 && fp <= fp_old && fp_old < 0f64 {
            break;
        }

        /* Check for maximum number of iterations */

        if iter == 10 {
            break;
        }

        /* Compute the Newton correction */

        compute_newton_correction(r, sdiag, perm, x, dxnorm, diag, w);

        /*#ifdef DEBUG
        printf ("newton_correction = ");
        gsl_vector_fprintf(stdout, w, "%g"); printf("\n");
        #endif*/

        {
            let wnorm = enorm(w);

            par_c = fp / (delta * wnorm * wnorm);
        }

        /*#ifdef DEBUG
        printf("fp = %g\n", fp);
        printf("par_lower = %g\n", par_lower);
        printf("par_upper = %g\n", par_upper);
        printf("par_c = %g\n", par_c);
        #endif*/


        /* Depending on the sign of the function, update par_lower or par_upper */

        if fp > 0f64 {
            if par > par_lower {
                par_lower = par;
                /*#ifdef DEBUG
                  printf("fp > 0: set par_lower = par = %g\n", par);
                #endif*/
            }
        } else if fp < 0f64 {
            if par < par_upper {
                /*#ifdef DEBUG
                  printf("fp < 0: set par_upper = par = %g\n", par);
                #endif*/
                par_upper = par;
            }
        }

        /* Compute an improved estimate for par */

        /*#ifdef DEBUG
          printf("improved estimate par = MAX(%g, %g) \n", par_lower, par+par_c);
        #endif*/

        par = par_lower.max(par + par_c);

        /*#ifdef DEBUG
          printf("improved estimate par = %g \n", par);
        #endif*/

    }

    /*#ifdef DEBUG
    printf("LEAVING lmpar, par = %g\n", par);
    #endif*/

    *par_inout = par;

    ::Value::Success
}

fn compute_newton_direction(r: &::MatrixF64, perm: &::Permutation, qtf: &::VectorF64, x: &mut ::VectorF64) {
    /* Compute and store in x the Gauss-Newton direction. If the
     Jacobian is rank-deficient then obtain a least squares
     solution. */

    let n = r.size2();

    for i in 0..n {
        let qtfi = qtf.get(i);

        x.set(i, qtfi);
    }

    let nsing = count_nsing(r);

    /*#ifdef DEBUG
    printf("nsing = %d\n", nsing);
    printf("r = "); gsl_matrix_fprintf(stdout, r, "%g"); printf("\n");
    printf("qtf = "); gsl_vector_fprintf(stdout, x, "%g"); printf("\n");
    #endif*/

    if nsing < n {
        for i in nsing..n {
            x.set(i, 0f64);
        }
    }

    if nsing > 0u64 {
        for j in nsing..0 {
            let rjj = r.get(j, j);
            let temp = x.get(j) / rjj;
          
            x.set(j, temp);
              
            for i in 0..j {
                let rij = r.get(i, j);
                let xi = x.get(i);
              
                x.set(i, xi - rij * temp);
            }
        }
    }

    perm.permute_vector_inverse(x);
}

fn compute_newton_bound(r: &::MatrixF64, x: &::VectorF64, dxnorm: f64, perm: &::Permutation, diag: &::VectorF64,
    w: &mut ::VectorF64) {
    /* If the jacobian is not rank-deficient then the Newton step
     provides a lower bound for the zero of the function. Otherwise
     set this bound to zero. */

    let n = r.size2();

    let nsing = count_nsing(r);

    if nsing < n {
        w.set_zero();
        return;
    }

    for i in 0..n {
        let pi = perm.get(i);

        let dpi = diag.get(pi);
        let xpi = x.get(pi);

        w.set(i, dpi * (dpi * xpi / dxnorm));
    }

    for j in 0..n {
        let mut sum = 0f64;

        for i in 0..j {
            sum += r.get(i, j) * w.get(i);
        }

        {
            let rjj = r.get(j, j);
            let wj = w.get(j);

            w.set(j, (wj - sum) / rjj);
        }
    }
}

/* This function computes the solution to the least squares system
   phi = [ A x =  b , lambda D x = 0 ]^2
    
   where A is an M by N matrix, D is an N by N diagonal matrix, lambda
   is a scalar parameter and b is a vector of length M.
   The function requires the factorization of A into A = Q R P^T,
   where Q is an orthogonal matrix, R is an upper triangular matrix
   with diagonal elements of non-increasing magnitude and P is a
   permuation matrix. The system above is then equivalent to
   [ R z = Q^T b, P^T (lambda D) P z = 0 ]
   where x = P z. If this system does not have full rank then a least
   squares solution is obtained.  On output the function also provides
   an upper triangular matrix S such that
   P^T (A^T A + lambda^2 D^T D) P = S^T S
   Parameters,
   
   r: On input, contains the full upper triangle of R. On output the
   strict lower triangle contains the transpose of the strict upper
   triangle of S, and the diagonal of S is stored in sdiag.  The full
   upper triangle of R is not modified.
   p: the encoded form of the permutation matrix P. column j of P is
   column p[j] of the identity matrix.
   lambda, diag: contains the scalar lambda and the diagonal elements
   of the matrix D
   qtb: contains the product Q^T b
   x: on output contains the least squares solution of the system
   wa: is a workspace of length N
   */
fn qrsolv(r: &mut ::MatrixF64, p: &::Permutation, lambda: f64, diag: &::VectorF64, qtb: &::VectorF64, 
    x: &mut ::VectorF64, sdiag: &mut ::VectorF64, wa: &mut ::VectorF64) -> ::Value {
    let n = r.size2();

    /* Copy r and qtb to preserve input and initialise s. In particular,
     save the diagonal elements of r in x */

    let mut j = 0;
    while j < n {
        let rjj = r.get(j, j);
        let qtbj = qtb.get(j);

        let mut i = j + 1;
        while i < n {
            let rji = r.get(j, i);

            r.set(i, j, rji);
            i += 1;
        }

        x.set(j, rjj);
        wa.set(j, qtbj);
        j += 1;
    }

    /* Eliminate the diagonal matrix d using a Givens rotation */
    j = 0;
    while j < n {
        let pj = p.get(j);

        let diagpj = lambda * diag.get(pj);

        if diagpj == 0f64 {
            continue;
        }

        sdiag.set(j, diagpj);

        let mut k = j + 1;
        while k < n {
            sdiag.set(k, 0f64);
            k += 1;
        }

        /* The transformations to eliminate the row of d modify only a
           single element of qtb beyond the first n, which is initially
           zero */

        let mut qtbpj = 0f64;
        k = j;
        while k < n {
            /* Determine a Givens rotation which eliminates the
               appropriate element in the current row of d */

            let wak = wa.get(k);
            let rkk = r.get(k, k);
            let sdiagk = sdiag.get(k);

            if sdiagk == 0f64 {
                continue;
            }

            let (sine, cosine) = if rkk.abs() < sdiagk.abs() {
                let cotangent = rkk / sdiagk;
                let t_sine = 0.5f64 / (0.25f64 + 0.25f64 * cotangent * cotangent).sqrt();

                (t_sine, t_sine * cotangent)
            } else {
                let tangent = sdiagk / rkk;
                let t_cos = 0.5f64 / (0.25f64 + 0.25f64 * tangent * tangent).sqrt();

                (t_cos * tangent, t_cos)
            };

            /* Compute the modified diagonal element of r and the
               modified element of [qtb,0] */

            {
                let new_rkk = cosine * rkk + sine * sdiagk;
                let new_wak = cosine * wak + sine * qtbpj;
                
                qtbpj = -sine * wak + cosine * qtbpj;

                r.set(k, k, new_rkk);
                wa.set(k, new_wak);
            }

            /* Accumulate the transformation in the row of s */
            let mut i = k + 1;
            while i < n {
                let rik = r.get(i, k);
                let sdiagi = sdiag.get(i);
              
                let new_rik = cosine * rik + sine * sdiagi;
                let new_sdiagi = -sine * rik + cosine * sdiagi;
              
                r.set(i, k, new_rik);
                sdiag.set(i, new_sdiagi);
                i += 1;
            }
            k += 1;
        }

        /* Store the corresponding diagonal element of s and restore the
           corresponding diagonal element of r */

        {
            let rjj = r.get(j, j);
            let xj = x.get(j);
            
            sdiag.set(j, rjj);
            r.set(j, j, xj);
        }
        j += 1;

    }

    /* Solve the triangular system for z. If the system is singular then
       obtain a least squares solution */

    let mut nsing = n;
    j = 0;
    while j < n {
        let sdiagj = sdiag.get(j);

        if sdiagj == 0f64 {
            nsing = j;
            break;
        }
        j += 1;
    }

    j = nsing;
    while j < n {
        wa.set(j, 0f64);
        j += 1;
    }

    let mut k = 0;
    while k < nsing {
        let mut sum = 0f64;

        j = (nsing - 1) - k;

        let mut i = j + 1;
        while i < nsing {
            sum += r.get(i, j) * wa.get(i);
            i += 1;
        }

        {
            let waj = wa.get(j);
            let sdiagj = sdiag.get(j);

            wa.set(j, (waj - sum) / sdiagj);
        }
        k += 1;
    }

    /* Permute the components of z back to the components of x */
    j = 0;
    while j < n {
        let pj = p.get(j);
        let waj = wa.get(j);

        x.set(pj, waj);
        j += 1;
    }

    ::Value::Success
}

fn compute_newton_correction(r: &::MatrixF64, sdiag: &::VectorF64, p: &::Permutation, x: &mut ::VectorF64, dxnorm: f64,
    diag: &::VectorF64, w: &mut ::VectorF64) {
    let n = r.size2();

    let mut i = 0;
    while i < n {
      let pi = p.get(i);

      let dpi = diag.get(pi);
      let xpi = x.get(pi);

      w.set(i, dpi * (dpi * xpi) / dxnorm);
      i += 1;
    }

    let mut j = 0;
    while j < n {
        let sj = sdiag.get(j);
        let wj = w.get(j);

        let tj = wj / sj;

        w.set(j, tj);

        let mut i = j + 1;
        while i < n {
            let rij = r.get (i, j);
            let wi = w.get (i);

            w.set (i, wi - rij * tj);
            i += 1;
        }
        j += 1;
    }
}

fn count_nsing(r: &::MatrixF64) -> u64 {
    /* Count the number of nonsingular entries. Returns the index of the
       first entry which is singular. */

    let n = r.size2();
    let mut i = 0u64;

    while i < n {
        let rii = r.get(i, i);

        if rii == 0f64 {
            break;
        }
        i += 1;
    }

    i
}