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

/*!
#Real Symmetric Matrices

For real symmetric matrices, the library uses the symmetric bidiagonalization and QR reduction
method. This is described in Golub & van Loan, section 8.3. The computed eigenvalues are accurate to
an absolute accuracy of \epsilon ||A||_2, where \epsilon is the machine precision.

#Complex Hermitian Matrices

For hermitian matrices, the library uses the complex form of the symmetric bidiagonalization and QR
reduction method.

#Real Nonsymmetric Matrices

The solution of the real nonsymmetric eigensystem problem for a matrix A involves computing the
Schur decomposition

A = Z T Z^T

where Z is an orthogonal matrix of Schur vectors and T, the Schur form, is quasi upper triangular
with diagonal 1-by-1 blocks which are real eigenvalues of A, and diagonal 2-by-2 blocks whose
eigenvalues are complex conjugate eigenvalues of A. The algorithm used is the double-shift Francis
method.

#Real Generalized Symmetric-Definite Eigensystems

The real generalized symmetric-definite eigenvalue problem is to find eigenvalues \lambda and
eigenvectors x such that

A x = lambda B x

where A and B are symmetric matrices, and B is positive-definite. This problem reduces to the
standard symmetric eigenvalue problem by applying the Cholesky decomposition to B:

```latex
                      A x = lambda B x
                      A x = lambda L L^t x
( L^{-1} A L^{-t} ) L^t x = lambda L^t x
```

Therefore, the problem becomes C y = lambda y where C = L^{-1} A L^{-t} is symmetric, and y = L^t x.
The standard symmetric eigensolver can be applied to the matrix C. The resulting eigenvectors are
backtransformed to find the vectors of the original problem. The eigenvalues and eigenvectors of the
generalized symmetric-definite eigenproblem are always real.

#Complex Generalized Hermitian-Definite Eigensystems

The complex generalized hermitian-definite eigenvalue problem is to find eigenvalues \lambda and
eigenvectors x such that

A x = \lambda B x

where A and B are hermitian matrices, and B is positive-definite. Similarly to the real case, this
can be reduced to C y = \lambda y where C = L^{-1} A L^{-H} is hermitian, and y = L^H x. The
standard hermitian eigensolver can be applied to the matrix C. The resulting eigenvectors are
backtransformed to find the vectors of the original problem. The eigenvalues of the generalized
hermitian-definite eigenproblem are always real.

#Real Generalized Nonsymmetric Eigensystems

Given two square matrices (A, B), the generalized nonsymmetric eigenvalue problem is to find
eigenvalues \lambda and eigenvectors x such that

A x = \lambda B x

We may also define the problem as finding eigenvalues \mu and eigenvectors y such that

\mu A y = B y
Note that these two problems are equivalent (with \lambda = 1/\mu) if neither \lambda nor \mu is
zero. If say, \lambda is zero, then it is still a well defined eigenproblem, but its alternate
problem involving \mu is not. Therefore, to allow for zero (and infinite) eigenvalues, the problem
which is actually solved is

\beta A x = \alpha B x
The eigensolver routines below will return two values \alpha and \beta and leave it to the user to
perform the divisions \lambda = \alpha / \beta and \mu = \beta / \alpha.

If the determinant of the matrix pencil A - \lambda B is zero for all \lambda, the problem is said
to be singular; otherwise it is called regular. Singularity normally leads to some
\alpha = \beta = 0 which means the eigenproblem is ill-conditioned and generally does not have well
defined eigenvalue solutions. The routines below are intended for regular matrix pencils and could
yield unpredictable results when applied to singular pencils.

The solution of the real generalized nonsymmetric eigensystem problem for a matrix pair (A, B)
involves computing the generalized Schur decomposition

A = Q S Z^T
B = Q T Z^T
where Q and Z are orthogonal matrices of left and right Schur vectors respectively, and (S, T) is
the generalized Schur form whose diagonal elements give the \alpha and \beta values. The algorithm
used is the QZ method due to Moler and Stewart (see references).
!*/

use ffi;
use enums;
use types::{MatrixF64, MatrixComplexF64, VectorF64, VectorComplexF64};

pub struct EigenSymmetricWorkspace {
    w: *mut ffi::gsl_eigen_symm_workspace,
}

impl EigenSymmetricWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n real symmetric
    /// matrices. The size of the workspace is O(2n).
    pub fn new(n: usize) -> Option<EigenSymmetricWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_symm_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenSymmetricWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues of the real symmetric matrix `A`. The diagonal and
    /// lower triangular part of `A` are destroyed during the computation, but the strict upper
    /// triangular part is not referenced. The eigenvalues are stored in the vector `eval` and are
    /// unordered.
    pub fn symm(&self, A: &mut MatrixF64, eval: &mut VectorF64) -> enums::Value {
        unsafe {
            ffi::gsl_eigen_symm(ffi::FFI::unwrap_unique(A),
                                ffi::FFI::unwrap_unique(eval),
                                self.w)
        }
    }
}

impl Drop for EigenSymmetricWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_symm_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_symm_workspace> for EigenSymmetricWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_symm_workspace) -> EigenSymmetricWorkspace {
        EigenSymmetricWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_symm_workspace) -> EigenSymmetricWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenSymmetricWorkspace) -> *const ffi::gsl_eigen_symm_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenSymmetricWorkspace) -> *mut ffi::gsl_eigen_symm_workspace {
        t.w
    }
}

pub struct EigenSymmetricVWorkspace {
    w: *mut ffi::gsl_eigen_symmv_workspace,
}

impl EigenSymmetricVWorkspace {
    /// This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n
    /// real symmetric matrices. The size of the workspace is O(4n).
    pub fn new(n: usize) -> Option<EigenSymmetricVWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_symmv_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenSymmetricVWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues and eigenvectors of the real symmetric matrix `A`.
    /// The diagonal and lower triangular part of `A` are destroyed during the computation, but the
    /// strict upper triangular part is not referenced. The eigenvalues are stored in the vector
    /// `eval` and are unordered. The corresponding eigenvectors are stored in the columns of the
    /// matrix `evec`. For example, the eigenvector in the first column corresponds to the first
    /// eigenvalue. The eigenvectors are guaranteed to be mutually orthogonal and normalised to unit
    /// magnitude.
    pub fn symmv(&self,
                 A: &mut MatrixF64,
                 eval: &mut VectorF64,
                 evec: &mut MatrixF64)
                 -> enums::Value {
        unsafe {
            ffi::gsl_eigen_symmv(ffi::FFI::unwrap_unique(A),
                                 ffi::FFI::unwrap_unique(eval),
                                 ffi::FFI::unwrap_unique(evec),
                                 self.w)
        }
    }
}

impl Drop for EigenSymmetricVWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_symmv_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_symmv_workspace> for EigenSymmetricVWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_symmv_workspace) -> EigenSymmetricVWorkspace {
        EigenSymmetricVWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_symmv_workspace) -> EigenSymmetricVWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenSymmetricVWorkspace) -> *const ffi::gsl_eigen_symmv_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenSymmetricVWorkspace) -> *mut ffi::gsl_eigen_symmv_workspace {
        t.w
    }
}

pub struct EigenHermitianWorkspace {
    w: *mut ffi::gsl_eigen_herm_workspace,
}

impl EigenHermitianWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n complex hermitian
    /// matrices. The size of the workspace is O(3n).
    pub fn new(n: usize) -> Option<EigenHermitianWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_herm_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenHermitianWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues of the complex hermitian matrix `A`. Additional
    /// workspace of the appropriate size must be provided in `self`. The diagonal and lower
    /// triangular part of `A` are destroyed during the computation, but the strict upper triangular
    /// part is not referenced. The imaginary parts of the diagonal are assumed to be zero and are
    /// not referenced. The eigenvalues are stored in the vector `eval` and are unordered.
    pub fn herm(&mut self, A: &mut MatrixComplexF64, eval: &mut VectorF64) -> enums::Value {
        unsafe {
            ffi::gsl_eigen_herm(ffi::FFI::unwrap_unique(A),
                                ffi::FFI::unwrap_unique(eval),
                                self.w)
        }
    }
}

impl Drop for EigenHermitianWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_herm_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_herm_workspace> for EigenHermitianWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_herm_workspace) -> EigenHermitianWorkspace {
        EigenHermitianWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_herm_workspace) -> EigenHermitianWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenHermitianWorkspace) -> *const ffi::gsl_eigen_herm_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenHermitianWorkspace) -> *mut ffi::gsl_eigen_herm_workspace {
        t.w
    }
}

pub struct EigenHermitianVWorkspace {
    w: *mut ffi::gsl_eigen_hermv_workspace,
}

impl EigenHermitianVWorkspace {
    /// This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n
    /// complex hermitian matrices. The size of the workspace is O(5n).
    pub fn new(n: usize) -> Option<EigenHermitianVWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_hermv_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenHermitianVWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues and eigenvectors of the complex hermitian matrix `A`.
    /// Additional workspace of the appropriate size must be provided in `self`. The diagonal and
    /// lower triangular part of `A` are destroyed during the computation, but the strict upper
    /// triangular part is not referenced. The imaginary parts of the diagonal are assumed to be
    /// zero and are not referenced. The eigenvalues are stored in the vector `eval` and are
    /// unordered. The corresponding complex eigenvectors are stored in the columns of the matrix
    /// `evec`. For example, the eigenvector in the first column corresponds to the first
    /// eigenvalue. The eigenvectors are guaranteed to be mutually orthogonal and normalised to unit
    /// magnitude.
    pub fn hermv(&mut self,
                 A: &mut MatrixComplexF64,
                 eval: &mut VectorF64,
                 evec: &mut MatrixComplexF64)
                 -> enums::Value {
        unsafe {
            ffi::gsl_eigen_hermv(ffi::FFI::unwrap_unique(A),
                                 ffi::FFI::unwrap_unique(eval),
                                 ffi::FFI::unwrap_unique(evec),
                                 self.w)
        }
    }
}

impl Drop for EigenHermitianVWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_hermv_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_hermv_workspace> for EigenHermitianVWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_hermv_workspace) -> EigenHermitianVWorkspace {
        EigenHermitianVWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_hermv_workspace) -> EigenHermitianVWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenHermitianVWorkspace) -> *const ffi::gsl_eigen_hermv_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenHermitianVWorkspace) -> *mut ffi::gsl_eigen_hermv_workspace {
        t.w
    }
}

pub struct EigenNonSymmWorkspace {
    w: *mut ffi::gsl_eigen_nonsymm_workspace,
}

impl EigenNonSymmWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n complex hermitian
    /// matrices. The size of the workspace is O(3n).
    pub fn new(n: usize) -> Option<EigenNonSymmWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_nonsymm_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenNonSymmWorkspace { w: tmp })
        }
    }

    /// This function sets some parameters which determine how the eigenvalue problem is solved in
    /// subsequent calls to gsl_eigen_nonsymm.
    ///
    /// If compute_t is set to 1, the full Schur form T will be computed by gsl_eigen_nonsymm. If it
    /// is set to 0, T will not be computed (this is the default setting). Computing the full Schur
    /// form T requires approximately 1.5–2 times the number of flops.
    ///
    /// If balance is set to 1, a balancing transformation is applied to the matrix prior to
    /// computing eigenvalues. This transformation is designed to make the rows and columns of the
    /// matrix have comparable norms, and can result in more accurate eigenvalues for matrices whose
    /// entries vary widely in magnitude. See
    /// [`Balancing`](http://www.gnu.org/software/gsl/manual/html_node/Balancing.html#Balancing) for
    /// more information. Note that the balancing transformation does not preserve the orthogonality
    /// of the Schur vectors, so if you wish to compute the Schur vectors with gsl_eigen_nonsymm_Z
    /// you will obtain the Schur vectors of the balanced matrix instead of the original matrix.
    /// The relationship will be
    ///
    /// T = Q^t D^(-1) A D Q
    ///
    /// where Q is the matrix of Schur vectors for the balanced matrix, and D is the balancing
    /// transformation. Then gsl_eigen_nonsymm_Z will compute a matrix Z which satisfies
    ///
    /// T = Z^(-1) A Z
    ///
    /// with Z = D Q. Note that Z will not be orthogonal. For this reason, balancing is not
    /// performed by default.
    pub fn params(&mut self, compute_t: i32, balance: i32) {
        unsafe { ffi::gsl_eigen_nonsymm_params(compute_t, balance, self.w) }
    }

    /// This function computes the eigenvalues of the real nonsymmetric matrix `A` and stores them
    /// in the vector `eval`. If T is desired, it is stored in the upper portion of `A` on output.
    /// Otherwise, on output, the diagonal of `A` will contain the 1-by-1 real eigenvalues and
    /// 2-by-2 complex conjugate eigenvalue systems, and the rest of `A` is destroyed. In rare
    /// cases, this function may fail to find all eigenvalues. If this happens, an error code is
    /// returned and the number of converged eigenvalues is stored in w->n_evals. The converged
    /// eigenvalues are stored in the beginning of `eval`.
    pub fn nonsymm(&mut self, A: &mut MatrixF64, eval: &mut VectorComplexF64) -> enums::Value {
        unsafe {
            ffi::gsl_eigen_nonsymm(ffi::FFI::unwrap_unique(A),
                                   ffi::FFI::unwrap_unique(eval),
                                   self.w)
        }
    }

    /// This function is identical to gsl_eigen_nonsymm except that it also computes the Schur
    /// vectors and stores them into `Z`.
    pub fn nonsymm_Z(&mut self,
                     A: &mut MatrixF64,
                     eval: &mut VectorComplexF64,
                     Z: &mut MatrixF64)
                     -> enums::Value {
        unsafe {
            ffi::gsl_eigen_nonsymm_Z(ffi::FFI::unwrap_unique(A),
                                     ffi::FFI::unwrap_unique(eval),
                                     ffi::FFI::unwrap_unique(Z),
                                     self.w)
        }
    }

    pub fn n_evals(&self) -> usize {
        unsafe { (*self.w).n_evals }
    }
}

impl Drop for EigenNonSymmWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_nonsymm_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_nonsymm_workspace> for EigenNonSymmWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_nonsymm_workspace) -> EigenNonSymmWorkspace {
        EigenNonSymmWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_nonsymm_workspace) -> EigenNonSymmWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenNonSymmWorkspace) -> *const ffi::gsl_eigen_nonsymm_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenNonSymmWorkspace) -> *mut ffi::gsl_eigen_nonsymm_workspace {
        t.w
    }
}

pub struct EigenNonSymmVWorkspace {
    w: *mut ffi::gsl_eigen_nonsymmv_workspace,
}

impl EigenNonSymmVWorkspace {
    /// This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n
    /// real nonsymmetric matrices. The size of the workspace is O(5n).
    pub fn new(n: usize) -> Option<EigenNonSymmVWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_nonsymmv_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenNonSymmVWorkspace { w: tmp })
        }
    }

    /// This function sets parameters which determine how the eigenvalue problem is solved in
    /// subsequent calls to gsl_eigen_nonsymmv. If `balance` is set to 1, a balancing transformation
    /// is applied to the matrix. See gsl_eigen_nonsymm_params for more information. Balancing is
    /// turned off by default since it does not preserve the orthogonality of the Schur vectors.
    pub fn params(&mut self, balance: i32) {
        unsafe { ffi::gsl_eigen_nonsymmv_params(balance, self.w) }
    }

    /// This function computes eigenvalues and right eigenvectors of the n-by-n real nonsymmetric
    /// matrix `A`. It first calls gsl_eigen_nonsymm to compute the eigenvalues, Schur form T, and
    /// Schur vectors. Then it finds eigenvectors of T and backtransforms them using the Schur
    /// vectors. The Schur vectors are destroyed in the process, but can be saved by using
    /// gsl_eigen_nonsymmv_Z. The computed eigenvectors are normalized to have unit magnitude. On
    /// output, the upper portion of `A` contains the Schur form T. If gsl_eigen_nonsymm fails, no
    /// eigenvectors are computed, and an error code is returned.
    pub fn nonsymmv(&mut self,
                    A: &mut MatrixF64,
                    eval: &mut VectorComplexF64,
                    evec: &mut MatrixComplexF64)
                    -> enums::Value {
        unsafe {
            ffi::gsl_eigen_nonsymmv(ffi::FFI::unwrap_unique(A),
                                    ffi::FFI::unwrap_unique(eval),
                                    ffi::FFI::unwrap_unique(evec),
                                    self.w)
        }
    }

    /// This function is identical to gsl_eigen_nonsymmv except that it also saves the Schur vectors
    /// into `Z`.
    pub fn nonsymmv_Z(&mut self,
                      A: &mut MatrixF64,
                      eval: &mut VectorComplexF64,
                      evec: &mut MatrixComplexF64,
                      Z: &mut MatrixF64)
                      -> enums::Value {
        unsafe {
            ffi::gsl_eigen_nonsymmv_Z(ffi::FFI::unwrap_unique(A),
                                      ffi::FFI::unwrap_unique(eval),
                                      ffi::FFI::unwrap_unique(evec),
                                      ffi::FFI::unwrap_unique(Z),
                                      self.w)
        }
    }
}

impl Drop for EigenNonSymmVWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_nonsymmv_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_nonsymmv_workspace> for EigenNonSymmVWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_nonsymmv_workspace) -> EigenNonSymmVWorkspace {
        EigenNonSymmVWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_nonsymmv_workspace) -> EigenNonSymmVWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenNonSymmVWorkspace) -> *const ffi::gsl_eigen_nonsymmv_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenNonSymmVWorkspace) -> *mut ffi::gsl_eigen_nonsymmv_workspace {
        t.w
    }
}

pub struct EigenGenSymmWorkspace {
    w: *mut ffi::gsl_eigen_gensymm_workspace,
}

impl EigenGenSymmWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n real generalized
    /// symmetric-definite eigensystems. The size of the workspace is O(2n).
    pub fn new(n: usize) -> Option<EigenGenSymmWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_gensymm_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenGenSymmWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues of the real generalized symmetric-definite matrix
    /// pair (A, B), and stores them in `eval`, using the method outlined above. On output, `B`
    /// contains its Cholesky decomposition and `A` is destroyed.
    pub fn gensymm(&mut self,
                   mut A: MatrixF64,
                   B: &mut MatrixF64,
                   eval: &mut VectorF64)
                   -> enums::Value {
        unsafe {
            ffi::gsl_eigen_gensymm(ffi::FFI::unwrap_unique(&mut A),
                                   ffi::FFI::unwrap_unique(B),
                                   ffi::FFI::unwrap_unique(eval),
                                   self.w)
        }
    }
}

impl Drop for EigenGenSymmWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_gensymm_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_gensymm_workspace> for EigenGenSymmWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_gensymm_workspace) -> EigenGenSymmWorkspace {
        EigenGenSymmWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_gensymm_workspace) -> EigenGenSymmWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenGenSymmWorkspace) -> *const ffi::gsl_eigen_gensymm_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenGenSymmWorkspace) -> *mut ffi::gsl_eigen_gensymm_workspace {
        t.w
    }
}

pub struct EigenGenSymmVWorkspace {
    w: *mut ffi::gsl_eigen_gensymmv_workspace,
}

impl EigenGenSymmVWorkspace {
    /// This function allocates a workspace for computing eigenvalues and eigenvectors of n-by-n
    /// real generalized symmetric-definite eigensystems. The size of the workspace is O(4n).
    pub fn new(n: usize) -> Option<EigenGenSymmVWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_gensymmv_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenGenSymmVWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues and eigenvectors of the real generalized
    /// symmetric-definite matrix pair (A, B), and stores them in `eval` and `evec` respectively.
    /// The computed eigenvectors are normalized to have unit magnitude. On output, `B` contains its
    /// Cholesky decomposition and `A` is destroyed.
    pub fn gensymmv(&mut self,
                    mut A: MatrixF64,
                    B: &mut MatrixF64,
                    eval: &mut VectorF64,
                    evec: &mut MatrixF64)
                    -> enums::Value {
        unsafe {
            ffi::gsl_eigen_gensymmv(ffi::FFI::unwrap_unique(&mut A),
                                    ffi::FFI::unwrap_unique(B),
                                    ffi::FFI::unwrap_unique(eval),
                                    ffi::FFI::unwrap_unique(evec),
                                    self.w)
        }
    }
}

impl Drop for EigenGenSymmVWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_gensymmv_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_gensymmv_workspace> for EigenGenSymmVWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_gensymmv_workspace) -> EigenGenSymmVWorkspace {
        EigenGenSymmVWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_gensymmv_workspace) -> EigenGenSymmVWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenGenSymmVWorkspace) -> *const ffi::gsl_eigen_gensymmv_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenGenSymmVWorkspace) -> *mut ffi::gsl_eigen_gensymmv_workspace {
        t.w
    }
}

pub struct EigenGenHermWorkspace {
    w: *mut ffi::gsl_eigen_genherm_workspace,
}

impl EigenGenHermWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n complex generalized
    /// hermitian-definite eigensystems. The size of the workspace is O(3n).
    pub fn new(n: usize) -> Option<EigenGenHermWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_genherm_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenGenHermWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues of the complex generalized hermitian-definite matrix
    /// pair (A, B), and stores them in `eval`, using the method outlined above. On output, `B`
    /// contains its Cholesky decomposition and `A` is destroyed.
    pub fn genherm(&mut self,
                   mut A: MatrixComplexF64,
                   B: &mut MatrixComplexF64,
                   eval: &mut VectorF64)
                   -> enums::Value {
        unsafe {
            ffi::gsl_eigen_genherm(ffi::FFI::unwrap_unique(&mut A),
                                   ffi::FFI::unwrap_unique(B),
                                   ffi::FFI::unwrap_unique(eval),
                                   self.w)
        }
    }
}

impl Drop for EigenGenHermWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_genherm_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_genherm_workspace> for EigenGenHermWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_genherm_workspace) -> EigenGenHermWorkspace {
        EigenGenHermWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_genherm_workspace) -> EigenGenHermWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenGenHermWorkspace) -> *const ffi::gsl_eigen_genherm_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenGenHermWorkspace) -> *mut ffi::gsl_eigen_genherm_workspace {
        t.w
    }
}

pub struct EigenGenHermVWorkspace {
    w: *mut ffi::gsl_eigen_genhermv_workspace,
}

impl EigenGenHermVWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n complex generalized
    /// hermitian-definite eigensystems. The size of the workspace is O(3n).
    pub fn new(n: usize) -> Option<EigenGenHermVWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_genhermv_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenGenHermVWorkspace { w: tmp })
        }
    }

    /// This function computes the eigenvalues of the complex generalized hermitian-definite matrix
    /// pair (A, B), and stores them in `eval`, using the method outlined above. On output, `B`
    /// contains its Cholesky decomposition and `A` is destroyed.
    pub fn genhermv(&mut self,
                    mut A: MatrixComplexF64,
                    B: &mut MatrixComplexF64,
                    eval: &mut VectorF64,
                    evec: &mut MatrixComplexF64)
                    -> enums::Value {
        unsafe {
            ffi::gsl_eigen_genhermv(ffi::FFI::unwrap_unique(&mut A),
                                    ffi::FFI::unwrap_unique(B),
                                    ffi::FFI::unwrap_unique(eval),
                                    ffi::FFI::unwrap_unique(evec),
                                    self.w)
        }
    }
}

impl Drop for EigenGenHermVWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_genhermv_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_genhermv_workspace> for EigenGenHermVWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_genhermv_workspace) -> EigenGenHermVWorkspace {
        EigenGenHermVWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_genhermv_workspace) -> EigenGenHermVWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenGenHermVWorkspace) -> *const ffi::gsl_eigen_genhermv_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenGenHermVWorkspace) -> *mut ffi::gsl_eigen_genhermv_workspace {
        t.w
    }
}

pub struct EigenGenWorkspace {
    w: *mut ffi::gsl_eigen_gen_workspace,
}

impl EigenGenWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n real generalized
    /// nonsymmetric eigensystems. The size of the workspace is O(n).
    pub fn new(n: usize) -> Option<EigenGenWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_gen_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenGenWorkspace { w: tmp })
        }
    }

    /// This function sets some parameters which determine how the eigenvalue problem is solved in
    /// subsequent calls to gsl_eigen_gen.
    ///
    /// If compute_s is set to 1, the full Schur form S will be computed by gsl_eigen_gen. If it is
    /// set to 0, S will not be computed (this is the default setting). S is a quasi upper
    /// triangular matrix with 1-by-1 and 2-by-2 blocks on its diagonal. 1-by-1 blocks correspond to
    /// real eigenvalues, and 2-by-2 blocks correspond to complex eigenvalues.
    ///
    /// If compute_t is set to 1, the full Schur form T will be computed by gsl_eigen_gen. If it is
    /// set to 0, T will not be computed (this is the default setting). T is an upper triangular
    /// matrix with non-negative elements on its diagonal. Any 2-by-2 blocks in S will correspond to
    /// a 2-by-2 diagonal block in T.
    ///
    /// The balance parameter is currently ignored, since generalized balancing is not yet
    /// implemented.
    pub fn params(&mut self, compute_s: i32, compute_t: i32, balance: i32) {
        unsafe { ffi::gsl_eigen_gen_params(compute_s, compute_t, balance, self.w) }
    }

    /// This function computes the eigenvalues of the real generalized nonsymmetric matrix pair
    /// (A, B), and stores them as pairs in (alpha, beta), where alpha is complex and beta is real.
    /// If \beta_i is non-zero, then \lambda = \alpha_i / \beta_i is an eigenvalue. Likewise, if
    /// \alpha_i is non-zero, then \mu = \beta_i / \alpha_i is an eigenvalue of the alternate
    /// problem \mu A y = B y. The elements of beta are normalized to be non-negative.
    ///
    /// If S is desired, it is stored in A on output. If T is desired, it is stored in B on output.
    /// The ordering of eigenvalues in (alpha, beta) follows the ordering of the diagonal blocks in
    /// the Schur forms S and T. In rare cases, this function may fail to find all eigenvalues. If
    /// this occurs, an error code is returned.
    pub fn gen(&mut self,
               A: &mut MatrixF64,
               B: &mut MatrixF64,
               alpha: &mut VectorComplexF64,
               beta: &mut VectorF64)
               -> enums::Value {
        unsafe {
            ffi::gsl_eigen_gen(ffi::FFI::unwrap_unique(A),
                               ffi::FFI::unwrap_unique(B),
                               ffi::FFI::unwrap_unique(alpha),
                               ffi::FFI::unwrap_unique(beta),
                               self.w)
        }
    }

    /// This function is identical to gsl_eigen_gen except that it also computes the left and right
    /// Schur vectors and stores them into `Q` and `Z` respectively.
    pub fn gen_QZ(&mut self,
                  A: &mut MatrixF64,
                  B: &mut MatrixF64,
                  alpha: &mut VectorComplexF64,
                  beta: &mut VectorF64,
                  Q: &mut MatrixF64,
                  Z: &mut MatrixF64)
                  -> enums::Value {
        unsafe {
            ffi::gsl_eigen_gen_QZ(ffi::FFI::unwrap_unique(A),
                                  ffi::FFI::unwrap_unique(B),
                                  ffi::FFI::unwrap_unique(alpha),
                                  ffi::FFI::unwrap_unique(beta),
                                  ffi::FFI::unwrap_unique(Q),
                                  ffi::FFI::unwrap_unique(Z),
                                  self.w)
        }
    }
}

impl Drop for EigenGenWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_gen_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_gen_workspace> for EigenGenWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_gen_workspace) -> EigenGenWorkspace {
        EigenGenWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_gen_workspace) -> EigenGenWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenGenWorkspace) -> *const ffi::gsl_eigen_gen_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenGenWorkspace) -> *mut ffi::gsl_eigen_gen_workspace {
        t.w
    }
}

pub struct EigenGenVWorkspace {
    w: *mut ffi::gsl_eigen_genv_workspace,
}

impl EigenGenVWorkspace {
    /// This function allocates a workspace for computing eigenvalues of n-by-n real generalized
    /// nonsymmetric eigensystems. The size of the workspace is O(n).
    pub fn new(n: usize) -> Option<EigenGenVWorkspace> {
        let tmp = unsafe { ffi::gsl_eigen_genv_alloc(n) };

        if tmp.is_null() {
            None
        } else {
            Some(EigenGenVWorkspace { w: tmp })
        }
    }

    /// This function computes eigenvalues and right eigenvectors of the n-by-n real generalized
    /// nonsymmetric matrix pair (A, B). The eigenvalues are stored in (alpha, beta) and the
    /// eigenvectors are stored in evec. It first calls gsl_eigen_gen to compute the eigenvalues,
    /// Schur forms, and Schur vectors. Then it finds eigenvectors of the Schur forms and
    /// backtransforms them using the Schur vectors. The Schur vectors are destroyed in the process,
    /// but can be saved by using gsl_eigen_genv_QZ. The computed eigenvectors are normalized to
    /// have unit magnitude. On output, (A, B) contains the generalized Schur form (S, T). If
    /// gsl_eigen_gen fails, no eigenvectors are computed, and an error code is returned.
    pub fn genv(&mut self,
                A: &mut MatrixF64,
                B: &mut MatrixF64,
                alpha: &mut VectorComplexF64,
                beta: &mut VectorF64,
                evec: &mut MatrixComplexF64)
                -> enums::Value {
        unsafe {
            ffi::gsl_eigen_genv(ffi::FFI::unwrap_unique(A),
                                ffi::FFI::unwrap_unique(B),
                                ffi::FFI::unwrap_unique(alpha),
                                ffi::FFI::unwrap_unique(beta),
                                ffi::FFI::unwrap_unique(evec),
                                self.w)
        }
    }

    /// This function is identical to gsl_eigen_genv except that it also computes the left and right
    /// Schur vectors and stores them into `Q` and `Z` respectively.
    pub fn genv_QZ(&mut self,
                   A: &mut MatrixF64,
                   B: &mut MatrixF64,
                   alpha: &mut VectorComplexF64,
                   beta: &mut VectorF64,
                   evec: &mut MatrixComplexF64,
                   Q: &mut MatrixF64,
                   Z: &mut MatrixF64)
                   -> enums::Value {
        unsafe {
            ffi::gsl_eigen_genv_QZ(ffi::FFI::unwrap_unique(A),
                                   ffi::FFI::unwrap_unique(B),
                                   ffi::FFI::unwrap_unique(alpha),
                                   ffi::FFI::unwrap_unique(beta),
                                   ffi::FFI::unwrap_unique(evec),
                                   ffi::FFI::unwrap_unique(Q),
                                   ffi::FFI::unwrap_unique(Z),
                                   self.w)
        }
    }
}

impl Drop for EigenGenVWorkspace {
    fn drop(&mut self) {
        unsafe { ffi::gsl_eigen_genv_free(self.w) };
        self.w = ::std::ptr::null_mut();
    }
}

impl ffi::FFI<ffi::gsl_eigen_genv_workspace> for EigenGenVWorkspace {
    fn wrap(t: *mut ffi::gsl_eigen_genv_workspace) -> EigenGenVWorkspace {
        EigenGenVWorkspace { w: t }
    }

    fn soft_wrap(t: *mut ffi::gsl_eigen_genv_workspace) -> EigenGenVWorkspace {
        Self::wrap(t)
    }

    fn unwrap_shared(t: &EigenGenVWorkspace) -> *const ffi::gsl_eigen_genv_workspace {
        t.w as *const _
    }

    fn unwrap_unique(t: &mut EigenGenVWorkspace) -> *mut ffi::gsl_eigen_genv_workspace {
        t.w
    }
}

#[test]
fn eigen_symmetric_workspace() {
    use MatrixF64;
    use VectorF64;

    let e = EigenSymmetricWorkspace::new(2).unwrap();
    let mut m = MatrixF64::new(2, 2).unwrap();

    let data = [5., 5., 1., 6.];
    m.set(0, 0, data[0]);
    m.set(0, 1, data[1]);
    m.set(1, 0, data[2]);
    m.set(1, 1, data[3]);
    let mut v = VectorF64::new(2).unwrap();
    e.symm(&mut m, &mut v);
    assert_eq!(&format!("{:.4} {:.4}", v.get(0), v.get(1)), "4.3820 6.6180");
}

// C code:
//
// ```ignore
// #include <gsl/gsl_eigen.h>
// #include <gsl/gsl_matrix.h>
//
// int main() {
//     gsl_eigen_symmv_workspace *t = gsl_eigen_symmv_alloc(3);
//     double data[] = {5., 5., 1., 6.};
//     gsl_matrix *m = gsl_matrix_calloc(2, 2);
//     gsl_matrix_set(m, 0, 0, data[0]);
//     gsl_matrix_set(m, 0, 1, data[1]);
//     gsl_matrix_set(m, 1, 0, data[2]);
//     gsl_matrix_set(m, 1, 1, data[3]);
//     gsl_vector *v = gsl_vector_calloc(2);
//     gsl_matrix *m2 = gsl_matrix_calloc(2, 2);
//     gsl_eigen_symmv(m, v, m2, t);
//     printf("%f %f\n", gsl_vector_get(v, 0), gsl_vector_get(v, 1));
//     printf("%f %f\n", gsl_matrix_get(m2, 0, 0), gsl_matrix_get(m2, 0, 1));
//     printf("%f %f\n", gsl_matrix_get(m2, 1, 0), gsl_matrix_get(m2, 1, 1));
//     gsl_eigen_symmv_free(t);
//     gsl_vector_free(v);
//     gsl_matrix_free(m);
//     gsl_matrix_free(m2);
//     return 0;
// }
// ```
#[test]
fn eigen_symmetric_vworkspace() {
    use MatrixF64;
    use VectorF64;

    let e = EigenSymmetricVWorkspace::new(3).unwrap();
    let data = [5., 5., 1., 6.];
    let mut m = MatrixF64::new(2, 2).unwrap();

    m.set(0, 0, data[0]);
    m.set(0, 1, data[1]);
    m.set(1, 0, data[2]);
    m.set(1, 1, data[3]);
    let mut m2 = MatrixF64::new(2, 2).unwrap();
    let mut v = VectorF64::new(2).unwrap();
    e.symmv(&mut m, &mut v, &mut m2);
    assert_eq!(&format!("{:.4} {:.4}", v.get(0), v.get(1)), "4.3820 6.6180");
    assert_eq!(&format!("{:.4} {:.4}", m2.get(0, 0), m2.get(0, 1)),
               "0.8507 0.5257");
    assert_eq!(&format!("{:.4} {:.4}", m2.get(1, 0), m2.get(1, 1)),
               "-0.5257 0.8507");
}

// C code:
//
// ```ignore
// #include <gsl/gsl_eigen.h>
// #include <gsl/gsl_matrix.h>
// #include <gsl/gsl_complex_math.h>
//
// int main() {
//     gsl_eigen_herm_workspace *t = gsl_eigen_herm_alloc(3);
//     gsl_matrix_complex *m = gsl_matrix_complex_calloc(2, 2);
//     gsl_complex c1 = gsl_complex_rect(5., 5.);
//     gsl_complex c2 = gsl_complex_rect(1., 4.);
//     gsl_complex c3 = gsl_complex_rect(2., 3.);
//     gsl_complex c4 = gsl_complex_rect(5., 7.);
//     gsl_matrix_complex_set(m, 0, 0, c1);
//     gsl_matrix_complex_set(m, 0, 1, c2);
//     gsl_matrix_complex_set(m, 1, 0, c3);
//     gsl_matrix_complex_set(m, 1, 1, c4);
//     gsl_vector *v = gsl_vector_calloc(2);
//     gsl_eigen_herm(m, v, t);
//     printf("%f %f\n", gsl_vector_get(v, 0), gsl_vector_get(v, 1));
//     gsl_eigen_herm_free(t);
//     gsl_vector_free(v);
//     gsl_matrix_complex_free(m);
//     return 0;
// }
// ```
#[test]
fn eigen_hermitian_workspace() {
    use MatrixComplexF64;
    use VectorF64;
    use ComplexF64;

    let mut e = EigenHermitianWorkspace::new(3).unwrap();
    let mut m = MatrixComplexF64::new(2, 2).unwrap();

    m.set(0, 0, &ComplexF64::rect(5., 5.));
    m.set(0, 1, &ComplexF64::rect(1., 4.));
    m.set(1, 0, &ComplexF64::rect(2., 3.));
    m.set(1, 1, &ComplexF64::rect(5., 7.));

    let mut v = VectorF64::new(2).unwrap();
    e.herm(&mut m, &mut v);
    assert_eq!(&format!("{:.4} {:.4}", v.get(0), v.get(1)), "8.6056 1.3944");
}

// C code:
//
// ```ignore
// #include <gsl/gsl_eigen.h>
// #include <gsl/gsl_matrix.h>
// #include <gsl/gsl_complex_math.h>
//
// int main() {
//     gsl_eigen_hermv_workspace *t = gsl_eigen_hermv_alloc(3);
//     gsl_matrix_complex *m = gsl_matrix_complex_calloc(2, 2);
//     gsl_complex c1 = gsl_complex_rect(5., 5.);
//     gsl_complex c2 = gsl_complex_rect(1., 4.);
//     gsl_complex c3 = gsl_complex_rect(2., 3.);
//     gsl_complex c4 = gsl_complex_rect(5., 7.);
//     gsl_matrix_complex_set(m, 0, 0, c1);
//     gsl_matrix_complex_set(m, 0, 1, c2);
//     gsl_matrix_complex_set(m, 1, 0, c3);
//     gsl_matrix_complex_set(m, 1, 1, c4);
//     gsl_vector *v = gsl_vector_calloc(2);
//     gsl_matrix_complex *m2 = gsl_matrix_complex_calloc(2, 2);
//     gsl_eigen_hermv(m, v, m2, t);
//     printf("%f %f\n", gsl_vector_get(v, 0), gsl_vector_get(v, 1));
//     printf("(%f, %f) (%f, %f)\n",
//            gsl_matrix_complex_get(m2, 0, 0).dat[0], gsl_matrix_complex_get(m2, 0, 0).dat[1],
//            gsl_matrix_complex_get(m2, 0, 1).dat[0], gsl_matrix_complex_get(m2, 0, 1).dat[1]);
//     printf("(%f, %f) (%f, %f)\n",
//            gsl_matrix_complex_get(m2, 1, 0).dat[0], gsl_matrix_complex_get(m2, 1, 0).dat[1],
//            gsl_matrix_complex_get(m2, 1, 1).dat[0], gsl_matrix_complex_get(m2, 1, 1).dat[1]);
//     gsl_eigen_hermv_free(t);
//     gsl_vector_free(v);
//     gsl_matrix_complex_free(m);
//     gsl_matrix_complex_free(m2);
//     return 0;
// }
// ```
#[test]
fn eigen_hermitian_vworkspace() {
    use ComplexF64;

    let mut e = EigenHermitianVWorkspace::new(3).unwrap();
    let mut m = MatrixComplexF64::new(2, 2).unwrap();

    m.set(0, 0, &ComplexF64::rect(5., 5.));
    m.set(0, 1, &ComplexF64::rect(1., 4.));
    m.set(1, 0, &ComplexF64::rect(2., 3.));
    m.set(1, 1, &ComplexF64::rect(5., 7.));

    let mut v = VectorF64::new(2).unwrap();
    let mut m2 = MatrixComplexF64::new(2, 2).unwrap();
    e.hermv(&mut m, &mut v, &mut m2);
    assert_eq!(&format!("{:.4} {:.4}", v.get(0), v.get(1)), "8.6056 1.3944");
    assert_eq!(&format!("({:.4}, {:.4}) ({:.4}, {:.4})",
                        m2.get(0, 0).dat[0],
                        m2.get(0, 0).dat[1],
                        m2.get(0, 1).dat[0],
                        m2.get(0, 1).dat[1]),
               "(0.7071, 0.0000) (0.7071, 0.0000)");
    assert_eq!(&format!("({:.4}, {:.4}) ({:.4}, {:.4})",
                        m2.get(1, 0).dat[0],
                        m2.get(1, 0).dat[1],
                        m2.get(1, 1).dat[0],
                        m2.get(1, 1).dat[1]),
               "(0.3922, 0.5883) (-0.3922, -0.5883)");
}