aboutsummaryrefslogtreecommitdiff
path: root/test/eigen2/gsl_helper.h
blob: d1d8545333a0634a84b76cce43ec43c1477445a5 (plain)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra. Eigen itself is part of the KDE project.
//
// Copyright (C) 2008 Gael Guennebaud <g.gael@free.fr>
//
// This Source Code Form is subject to the terms of the Mozilla
// Public License v. 2.0. If a copy of the MPL was not distributed
// with this file, You can obtain one at http://mozilla.org/MPL/2.0/.

#ifndef EIGEN_GSL_HELPER
#define EIGEN_GSL_HELPER

#include <Eigen/Core>

#include <gsl/gsl_blas.h>
#include <gsl/gsl_multifit.h>
#include <gsl/gsl_eigen.h>
#include <gsl/gsl_linalg.h>
#include <gsl/gsl_complex.h>
#include <gsl/gsl_complex_math.h>

namespace Eigen {

template<typename Scalar, bool IsComplex = NumTraits<Scalar>::IsComplex> struct GslTraits
{
  typedef gsl_matrix* Matrix;
  typedef gsl_vector* Vector;
  static Matrix createMatrix(int rows, int cols) { return gsl_matrix_alloc(rows,cols); }
  static Vector createVector(int size) { return gsl_vector_alloc(size); }
  static void free(Matrix& m) { gsl_matrix_free(m); m=0; }
  static void free(Vector& m) { gsl_vector_free(m); m=0; }
  static void prod(const Matrix& m, const Vector& v, Vector& x) { gsl_blas_dgemv(CblasNoTrans,1,m,v,0,x); }
  static void cholesky(Matrix& m) { gsl_linalg_cholesky_decomp(m); }
  static void cholesky_solve(const Matrix& m, const Vector& b, Vector& x) { gsl_linalg_cholesky_solve(m,b,x); }
  static void eigen_symm(const Matrix& m, Vector& eval, Matrix& evec)
  {
    gsl_eigen_symmv_workspace * w = gsl_eigen_symmv_alloc(m->size1);
    Matrix a = createMatrix(m->size1, m->size2);
    gsl_matrix_memcpy(a, m);
    gsl_eigen_symmv(a,eval,evec,w);
    gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
    gsl_eigen_symmv_free(w);
    free(a);
  }
  static void eigen_symm_gen(const Matrix& m, const Matrix& _b, Vector& eval, Matrix& evec)
  {
    gsl_eigen_gensymmv_workspace * w = gsl_eigen_gensymmv_alloc(m->size1);
    Matrix a = createMatrix(m->size1, m->size2);
    Matrix b = createMatrix(_b->size1, _b->size2);
    gsl_matrix_memcpy(a, m);
    gsl_matrix_memcpy(b, _b);
    gsl_eigen_gensymmv(a,b,eval,evec,w);
    gsl_eigen_symmv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
    gsl_eigen_gensymmv_free(w);
    free(a);
  }
};

template<typename Scalar> struct GslTraits<Scalar,true>
{
  typedef gsl_matrix_complex* Matrix;
  typedef gsl_vector_complex* Vector;
  static Matrix createMatrix(int rows, int cols) { return gsl_matrix_complex_alloc(rows,cols); }
  static Vector createVector(int size) { return gsl_vector_complex_alloc(size); }
  static void free(Matrix& m) { gsl_matrix_complex_free(m); m=0; }
  static void free(Vector& m) { gsl_vector_complex_free(m); m=0; }
  static void cholesky(Matrix& m) { gsl_linalg_complex_cholesky_decomp(m); }
  static void cholesky_solve(const Matrix& m, const Vector& b, Vector& x) { gsl_linalg_complex_cholesky_solve(m,b,x); }
  static void prod(const Matrix& m, const Vector& v, Vector& x)
  { gsl_blas_zgemv(CblasNoTrans,gsl_complex_rect(1,0),m,v,gsl_complex_rect(0,0),x); }
  static void eigen_symm(const Matrix& m, gsl_vector* &eval, Matrix& evec)
  {
    gsl_eigen_hermv_workspace * w = gsl_eigen_hermv_alloc(m->size1);
    Matrix a = createMatrix(m->size1, m->size2);
    gsl_matrix_complex_memcpy(a, m);
    gsl_eigen_hermv(a,eval,evec,w);
    gsl_eigen_hermv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
    gsl_eigen_hermv_free(w);
    free(a);
  }
  static void eigen_symm_gen(const Matrix& m, const Matrix& _b, gsl_vector* &eval, Matrix& evec)
  {
    gsl_eigen_genhermv_workspace * w = gsl_eigen_genhermv_alloc(m->size1);
    Matrix a = createMatrix(m->size1, m->size2);
    Matrix b = createMatrix(_b->size1, _b->size2);
    gsl_matrix_complex_memcpy(a, m);
    gsl_matrix_complex_memcpy(b, _b);
    gsl_eigen_genhermv(a,b,eval,evec,w);
    gsl_eigen_hermv_sort(eval, evec, GSL_EIGEN_SORT_VAL_ASC);
    gsl_eigen_genhermv_free(w);
    free(a);
  }
};

template<typename MatrixType>
void convert(const MatrixType& m, gsl_matrix* &res)
{
//   if (res)
//     gsl_matrix_free(res);
  res = gsl_matrix_alloc(m.rows(), m.cols());
  for (int i=0 ; i<m.rows() ; ++i)
    for (int j=0 ; j<m.cols(); ++j)
      gsl_matrix_set(res, i, j, m(i,j));
}

template<typename MatrixType>
void convert(const gsl_matrix* m, MatrixType& res)
{
  res.resize(int(m->size1), int(m->size2));
  for (int i=0 ; i<res.rows() ; ++i)
    for (int j=0 ; j<res.cols(); ++j)
      res(i,j) = gsl_matrix_get(m,i,j);
}

template<typename VectorType>
void convert(const VectorType& m, gsl_vector* &res)
{
  if (res) gsl_vector_free(res);
  res = gsl_vector_alloc(m.size());
  for (int i=0 ; i<m.size() ; ++i)
      gsl_vector_set(res, i, m[i]);
}

template<typename VectorType>
void convert(const gsl_vector* m, VectorType& res)
{
  res.resize (m->size);
  for (int i=0 ; i<res.rows() ; ++i)
    res[i] = gsl_vector_get(m, i);
}

template<typename MatrixType>
void convert(const MatrixType& m, gsl_matrix_complex* &res)
{
  res = gsl_matrix_complex_alloc(m.rows(), m.cols());
  for (int i=0 ; i<m.rows() ; ++i)
    for (int j=0 ; j<m.cols(); ++j)
    {
      gsl_matrix_complex_set(res, i, j,
        gsl_complex_rect(m(i,j).real(), m(i,j).imag()));
    }
}

template<typename MatrixType>
void convert(const gsl_matrix_complex* m, MatrixType& res)
{
  res.resize(int(m->size1), int(m->size2));
  for (int i=0 ; i<res.rows() ; ++i)
    for (int j=0 ; j<res.cols(); ++j)
      res(i,j) = typename MatrixType::Scalar(
        GSL_REAL(gsl_matrix_complex_get(m,i,j)),
        GSL_IMAG(gsl_matrix_complex_get(m,i,j)));
}

template<typename VectorType>
void convert(const VectorType& m, gsl_vector_complex* &res)
{
  res = gsl_vector_complex_alloc(m.size());
  for (int i=0 ; i<m.size() ; ++i)
      gsl_vector_complex_set(res, i, gsl_complex_rect(m[i].real(), m[i].imag()));
}

template<typename VectorType>
void convert(const gsl_vector_complex* m, VectorType& res)
{
  res.resize(m->size);
  for (int i=0 ; i<res.rows() ; ++i)
    res[i] = typename VectorType::Scalar(
        GSL_REAL(gsl_vector_complex_get(m, i)),
        GSL_IMAG(gsl_vector_complex_get(m, i)));
}

}

#endif // EIGEN_GSL_HELPER