aboutsummaryrefslogtreecommitdiffhomepage
path: root/demos/mix_eigen_and_c/binary_library.cpp
blob: 15a2d03e9ae8adeb20bff36269aa78fca1285520 (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
176
177
178
179
180
181
182
183
184
185
// This file is part of Eigen, a lightweight C++ template library
// for linear algebra.
//
// Copyright (C) 2009 Benoit Jacob <jacob.benoit.1@gmail.com>
//
// 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/.

// This C++ file compiles to binary code that can be linked to by your C program,
// thanks to the extern "C" syntax used in the declarations in binary_library.h.

#include "binary_library.h"

#include <Eigen/Core>

using namespace Eigen;

/************************* pointer conversion methods **********************************************/

////// class MatrixXd //////

inline MatrixXd& c_to_eigen(C_MatrixXd* ptr)
{
  return *reinterpret_cast<MatrixXd*>(ptr);
}

inline const MatrixXd& c_to_eigen(const C_MatrixXd* ptr)
{
  return *reinterpret_cast<const MatrixXd*>(ptr);
}

inline C_MatrixXd* eigen_to_c(MatrixXd& ref)
{
  return reinterpret_cast<C_MatrixXd*>(&ref);
}

inline const C_MatrixXd* eigen_to_c(const MatrixXd& ref)
{
  return reinterpret_cast<const C_MatrixXd*>(&ref);
}

////// class Map<MatrixXd> //////

inline Map<MatrixXd>& c_to_eigen(C_Map_MatrixXd* ptr)
{
  return *reinterpret_cast<Map<MatrixXd>*>(ptr);
}

inline const Map<MatrixXd>& c_to_eigen(const C_Map_MatrixXd* ptr)
{
  return *reinterpret_cast<const Map<MatrixXd>*>(ptr);
}

inline C_Map_MatrixXd* eigen_to_c(Map<MatrixXd>& ref)
{
  return reinterpret_cast<C_Map_MatrixXd*>(&ref);
}

inline const C_Map_MatrixXd* eigen_to_c(const Map<MatrixXd>& ref)
{
  return reinterpret_cast<const C_Map_MatrixXd*>(&ref);
}


/************************* implementation of classes **********************************************/


////// class MatrixXd //////


C_MatrixXd* MatrixXd_new(int rows, int cols)
{
  return eigen_to_c(*new MatrixXd(rows,cols));
}

void MatrixXd_delete(C_MatrixXd *m)
{
  delete &c_to_eigen(m);
}

double* MatrixXd_data(C_MatrixXd *m)
{
  return c_to_eigen(m).data();
}

void MatrixXd_set_zero(C_MatrixXd *m)
{
  c_to_eigen(m).setZero();
}

void MatrixXd_resize(C_MatrixXd *m, int rows, int cols)
{
  c_to_eigen(m).resize(rows,cols);
}

void MatrixXd_copy(C_MatrixXd *dst, const C_MatrixXd *src)
{
  c_to_eigen(dst) = c_to_eigen(src);
}

void MatrixXd_copy_map(C_MatrixXd *dst, const C_Map_MatrixXd *src)
{
  c_to_eigen(dst) = c_to_eigen(src);
}

void MatrixXd_set_coeff(C_MatrixXd *m, int i, int j, double coeff)
{
  c_to_eigen(m)(i,j) = coeff;
}

double MatrixXd_get_coeff(const C_MatrixXd *m, int i, int j)
{
  return c_to_eigen(m)(i,j);
}

void MatrixXd_print(const C_MatrixXd *m)
{
  std::cout << c_to_eigen(m) << std::endl;
}

void MatrixXd_multiply(const C_MatrixXd *m1, const C_MatrixXd *m2, C_MatrixXd *result)
{
  c_to_eigen(result) = c_to_eigen(m1) * c_to_eigen(m2);
}

void MatrixXd_add(const C_MatrixXd *m1, const C_MatrixXd *m2, C_MatrixXd *result)
{
  c_to_eigen(result) = c_to_eigen(m1) + c_to_eigen(m2);
}



////// class Map_MatrixXd //////


C_Map_MatrixXd* Map_MatrixXd_new(double *array, int rows, int cols)
{
  return eigen_to_c(*new Map<MatrixXd>(array,rows,cols));
}

void Map_MatrixXd_delete(C_Map_MatrixXd *m)
{
  delete &c_to_eigen(m);
}

void Map_MatrixXd_set_zero(C_Map_MatrixXd *m)
{
  c_to_eigen(m).setZero();
}

void Map_MatrixXd_copy(C_Map_MatrixXd *dst, const C_Map_MatrixXd *src)
{
  c_to_eigen(dst) = c_to_eigen(src);
}

void Map_MatrixXd_copy_matrix(C_Map_MatrixXd *dst, const C_MatrixXd *src)
{
  c_to_eigen(dst) = c_to_eigen(src);
}

void Map_MatrixXd_set_coeff(C_Map_MatrixXd *m, int i, int j, double coeff)
{
  c_to_eigen(m)(i,j) = coeff;
}

double Map_MatrixXd_get_coeff(const C_Map_MatrixXd *m, int i, int j)
{
  return c_to_eigen(m)(i,j);
}

void Map_MatrixXd_print(const C_Map_MatrixXd *m)
{
  std::cout << c_to_eigen(m) << std::endl;
}

void Map_MatrixXd_multiply(const C_Map_MatrixXd *m1, const C_Map_MatrixXd *m2, C_Map_MatrixXd *result)
{
  c_to_eigen(result) = c_to_eigen(m1) * c_to_eigen(m2);
}

void Map_MatrixXd_add(const C_Map_MatrixXd *m1, const C_Map_MatrixXd *m2, C_Map_MatrixXd *result)
{
  c_to_eigen(result) = c_to_eigen(m1) + c_to_eigen(m2);
}