1 // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9
10 // This C++ file compiles to binary code that can be linked to by your C program,
11 // thanks to the extern "C" syntax used in the declarations in binary_library.h.
12
13 #include "binary_library.h"
14
15 #include <Eigen/Core>
16
17 using namespace Eigen;
18
19 /************************* pointer conversion methods **********************************************/
20
21 ////// class MatrixXd //////
22
c_to_eigen(C_MatrixXd * ptr)23 inline MatrixXd& c_to_eigen(C_MatrixXd* ptr)
24 {
25 return *reinterpret_cast<MatrixXd*>(ptr);
26 }
27
c_to_eigen(const C_MatrixXd * ptr)28 inline const MatrixXd& c_to_eigen(const C_MatrixXd* ptr)
29 {
30 return *reinterpret_cast<const MatrixXd*>(ptr);
31 }
32
eigen_to_c(MatrixXd & ref)33 inline C_MatrixXd* eigen_to_c(MatrixXd& ref)
34 {
35 return reinterpret_cast<C_MatrixXd*>(&ref);
36 }
37
eigen_to_c(const MatrixXd & ref)38 inline const C_MatrixXd* eigen_to_c(const MatrixXd& ref)
39 {
40 return reinterpret_cast<const C_MatrixXd*>(&ref);
41 }
42
43 ////// class Map<MatrixXd> //////
44
c_to_eigen(C_Map_MatrixXd * ptr)45 inline Map<MatrixXd>& c_to_eigen(C_Map_MatrixXd* ptr)
46 {
47 return *reinterpret_cast<Map<MatrixXd>*>(ptr);
48 }
49
c_to_eigen(const C_Map_MatrixXd * ptr)50 inline const Map<MatrixXd>& c_to_eigen(const C_Map_MatrixXd* ptr)
51 {
52 return *reinterpret_cast<const Map<MatrixXd>*>(ptr);
53 }
54
eigen_to_c(Map<MatrixXd> & ref)55 inline C_Map_MatrixXd* eigen_to_c(Map<MatrixXd>& ref)
56 {
57 return reinterpret_cast<C_Map_MatrixXd*>(&ref);
58 }
59
eigen_to_c(const Map<MatrixXd> & ref)60 inline const C_Map_MatrixXd* eigen_to_c(const Map<MatrixXd>& ref)
61 {
62 return reinterpret_cast<const C_Map_MatrixXd*>(&ref);
63 }
64
65
66 /************************* implementation of classes **********************************************/
67
68
69 ////// class MatrixXd //////
70
71
MatrixXd_new(int rows,int cols)72 C_MatrixXd* MatrixXd_new(int rows, int cols)
73 {
74 return eigen_to_c(*new MatrixXd(rows,cols));
75 }
76
MatrixXd_delete(C_MatrixXd * m)77 void MatrixXd_delete(C_MatrixXd *m)
78 {
79 delete &c_to_eigen(m);
80 }
81
MatrixXd_data(C_MatrixXd * m)82 double* MatrixXd_data(C_MatrixXd *m)
83 {
84 return c_to_eigen(m).data();
85 }
86
MatrixXd_set_zero(C_MatrixXd * m)87 void MatrixXd_set_zero(C_MatrixXd *m)
88 {
89 c_to_eigen(m).setZero();
90 }
91
MatrixXd_resize(C_MatrixXd * m,int rows,int cols)92 void MatrixXd_resize(C_MatrixXd *m, int rows, int cols)
93 {
94 c_to_eigen(m).resize(rows,cols);
95 }
96
MatrixXd_copy(C_MatrixXd * dst,const C_MatrixXd * src)97 void MatrixXd_copy(C_MatrixXd *dst, const C_MatrixXd *src)
98 {
99 c_to_eigen(dst) = c_to_eigen(src);
100 }
101
MatrixXd_copy_map(C_MatrixXd * dst,const C_Map_MatrixXd * src)102 void MatrixXd_copy_map(C_MatrixXd *dst, const C_Map_MatrixXd *src)
103 {
104 c_to_eigen(dst) = c_to_eigen(src);
105 }
106
MatrixXd_set_coeff(C_MatrixXd * m,int i,int j,double coeff)107 void MatrixXd_set_coeff(C_MatrixXd *m, int i, int j, double coeff)
108 {
109 c_to_eigen(m)(i,j) = coeff;
110 }
111
MatrixXd_get_coeff(const C_MatrixXd * m,int i,int j)112 double MatrixXd_get_coeff(const C_MatrixXd *m, int i, int j)
113 {
114 return c_to_eigen(m)(i,j);
115 }
116
MatrixXd_print(const C_MatrixXd * m)117 void MatrixXd_print(const C_MatrixXd *m)
118 {
119 std::cout << c_to_eigen(m) << std::endl;
120 }
121
MatrixXd_multiply(const C_MatrixXd * m1,const C_MatrixXd * m2,C_MatrixXd * result)122 void MatrixXd_multiply(const C_MatrixXd *m1, const C_MatrixXd *m2, C_MatrixXd *result)
123 {
124 c_to_eigen(result) = c_to_eigen(m1) * c_to_eigen(m2);
125 }
126
MatrixXd_add(const C_MatrixXd * m1,const C_MatrixXd * m2,C_MatrixXd * result)127 void MatrixXd_add(const C_MatrixXd *m1, const C_MatrixXd *m2, C_MatrixXd *result)
128 {
129 c_to_eigen(result) = c_to_eigen(m1) + c_to_eigen(m2);
130 }
131
132
133
134 ////// class Map_MatrixXd //////
135
136
Map_MatrixXd_new(double * array,int rows,int cols)137 C_Map_MatrixXd* Map_MatrixXd_new(double *array, int rows, int cols)
138 {
139 return eigen_to_c(*new Map<MatrixXd>(array,rows,cols));
140 }
141
Map_MatrixXd_delete(C_Map_MatrixXd * m)142 void Map_MatrixXd_delete(C_Map_MatrixXd *m)
143 {
144 delete &c_to_eigen(m);
145 }
146
Map_MatrixXd_set_zero(C_Map_MatrixXd * m)147 void Map_MatrixXd_set_zero(C_Map_MatrixXd *m)
148 {
149 c_to_eigen(m).setZero();
150 }
151
Map_MatrixXd_copy(C_Map_MatrixXd * dst,const C_Map_MatrixXd * src)152 void Map_MatrixXd_copy(C_Map_MatrixXd *dst, const C_Map_MatrixXd *src)
153 {
154 c_to_eigen(dst) = c_to_eigen(src);
155 }
156
Map_MatrixXd_copy_matrix(C_Map_MatrixXd * dst,const C_MatrixXd * src)157 void Map_MatrixXd_copy_matrix(C_Map_MatrixXd *dst, const C_MatrixXd *src)
158 {
159 c_to_eigen(dst) = c_to_eigen(src);
160 }
161
Map_MatrixXd_set_coeff(C_Map_MatrixXd * m,int i,int j,double coeff)162 void Map_MatrixXd_set_coeff(C_Map_MatrixXd *m, int i, int j, double coeff)
163 {
164 c_to_eigen(m)(i,j) = coeff;
165 }
166
Map_MatrixXd_get_coeff(const C_Map_MatrixXd * m,int i,int j)167 double Map_MatrixXd_get_coeff(const C_Map_MatrixXd *m, int i, int j)
168 {
169 return c_to_eigen(m)(i,j);
170 }
171
Map_MatrixXd_print(const C_Map_MatrixXd * m)172 void Map_MatrixXd_print(const C_Map_MatrixXd *m)
173 {
174 std::cout << c_to_eigen(m) << std::endl;
175 }
176
Map_MatrixXd_multiply(const C_Map_MatrixXd * m1,const C_Map_MatrixXd * m2,C_Map_MatrixXd * result)177 void Map_MatrixXd_multiply(const C_Map_MatrixXd *m1, const C_Map_MatrixXd *m2, C_Map_MatrixXd *result)
178 {
179 c_to_eigen(result) = c_to_eigen(m1) * c_to_eigen(m2);
180 }
181
Map_MatrixXd_add(const C_Map_MatrixXd * m1,const C_Map_MatrixXd * m2,C_Map_MatrixXd * result)182 void Map_MatrixXd_add(const C_Map_MatrixXd *m1, const C_Map_MatrixXd *m2, C_Map_MatrixXd *result)
183 {
184 c_to_eigen(result) = c_to_eigen(m1) + c_to_eigen(m2);
185 }
186