binary_library.cpp 4.1 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185
  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. // This C++ file compiles to binary code that can be linked to by your C program,
  10. // thanks to the extern "C" syntax used in the declarations in binary_library.h.
  11. #include "binary_library.h"
  12. #include <Eigen/Core>
  13. using namespace Eigen;
  14. /************************* pointer conversion methods **********************************************/
  15. ////// class MatrixXd //////
  16. inline MatrixXd& c_to_eigen(C_MatrixXd* ptr)
  17. {
  18. return *reinterpret_cast<MatrixXd*>(ptr);
  19. }
  20. inline const MatrixXd& c_to_eigen(const C_MatrixXd* ptr)
  21. {
  22. return *reinterpret_cast<const MatrixXd*>(ptr);
  23. }
  24. inline C_MatrixXd* eigen_to_c(MatrixXd& ref)
  25. {
  26. return reinterpret_cast<C_MatrixXd*>(&ref);
  27. }
  28. inline const C_MatrixXd* eigen_to_c(const MatrixXd& ref)
  29. {
  30. return reinterpret_cast<const C_MatrixXd*>(&ref);
  31. }
  32. ////// class Map<MatrixXd> //////
  33. inline Map<MatrixXd>& c_to_eigen(C_Map_MatrixXd* ptr)
  34. {
  35. return *reinterpret_cast<Map<MatrixXd>*>(ptr);
  36. }
  37. inline const Map<MatrixXd>& c_to_eigen(const C_Map_MatrixXd* ptr)
  38. {
  39. return *reinterpret_cast<const Map<MatrixXd>*>(ptr);
  40. }
  41. inline C_Map_MatrixXd* eigen_to_c(Map<MatrixXd>& ref)
  42. {
  43. return reinterpret_cast<C_Map_MatrixXd*>(&ref);
  44. }
  45. inline const C_Map_MatrixXd* eigen_to_c(const Map<MatrixXd>& ref)
  46. {
  47. return reinterpret_cast<const C_Map_MatrixXd*>(&ref);
  48. }
  49. /************************* implementation of classes **********************************************/
  50. ////// class MatrixXd //////
  51. C_MatrixXd* MatrixXd_new(int rows, int cols)
  52. {
  53. return eigen_to_c(*new MatrixXd(rows,cols));
  54. }
  55. void MatrixXd_delete(C_MatrixXd *m)
  56. {
  57. delete &c_to_eigen(m);
  58. }
  59. double* MatrixXd_data(C_MatrixXd *m)
  60. {
  61. return c_to_eigen(m).data();
  62. }
  63. void MatrixXd_set_zero(C_MatrixXd *m)
  64. {
  65. c_to_eigen(m).setZero();
  66. }
  67. void MatrixXd_resize(C_MatrixXd *m, int rows, int cols)
  68. {
  69. c_to_eigen(m).resize(rows,cols);
  70. }
  71. void MatrixXd_copy(C_MatrixXd *dst, const C_MatrixXd *src)
  72. {
  73. c_to_eigen(dst) = c_to_eigen(src);
  74. }
  75. void MatrixXd_copy_map(C_MatrixXd *dst, const C_Map_MatrixXd *src)
  76. {
  77. c_to_eigen(dst) = c_to_eigen(src);
  78. }
  79. void MatrixXd_set_coeff(C_MatrixXd *m, int i, int j, double coeff)
  80. {
  81. c_to_eigen(m)(i,j) = coeff;
  82. }
  83. double MatrixXd_get_coeff(const C_MatrixXd *m, int i, int j)
  84. {
  85. return c_to_eigen(m)(i,j);
  86. }
  87. void MatrixXd_print(const C_MatrixXd *m)
  88. {
  89. std::cout << c_to_eigen(m) << std::endl;
  90. }
  91. void MatrixXd_multiply(const C_MatrixXd *m1, const C_MatrixXd *m2, C_MatrixXd *result)
  92. {
  93. c_to_eigen(result) = c_to_eigen(m1) * c_to_eigen(m2);
  94. }
  95. void MatrixXd_add(const C_MatrixXd *m1, const C_MatrixXd *m2, C_MatrixXd *result)
  96. {
  97. c_to_eigen(result) = c_to_eigen(m1) + c_to_eigen(m2);
  98. }
  99. ////// class Map_MatrixXd //////
  100. C_Map_MatrixXd* Map_MatrixXd_new(double *array, int rows, int cols)
  101. {
  102. return eigen_to_c(*new Map<MatrixXd>(array,rows,cols));
  103. }
  104. void Map_MatrixXd_delete(C_Map_MatrixXd *m)
  105. {
  106. delete &c_to_eigen(m);
  107. }
  108. void Map_MatrixXd_set_zero(C_Map_MatrixXd *m)
  109. {
  110. c_to_eigen(m).setZero();
  111. }
  112. void Map_MatrixXd_copy(C_Map_MatrixXd *dst, const C_Map_MatrixXd *src)
  113. {
  114. c_to_eigen(dst) = c_to_eigen(src);
  115. }
  116. void Map_MatrixXd_copy_matrix(C_Map_MatrixXd *dst, const C_MatrixXd *src)
  117. {
  118. c_to_eigen(dst) = c_to_eigen(src);
  119. }
  120. void Map_MatrixXd_set_coeff(C_Map_MatrixXd *m, int i, int j, double coeff)
  121. {
  122. c_to_eigen(m)(i,j) = coeff;
  123. }
  124. double Map_MatrixXd_get_coeff(const C_Map_MatrixXd *m, int i, int j)
  125. {
  126. return c_to_eigen(m)(i,j);
  127. }
  128. void Map_MatrixXd_print(const C_Map_MatrixXd *m)
  129. {
  130. std::cout << c_to_eigen(m) << std::endl;
  131. }
  132. void Map_MatrixXd_multiply(const C_Map_MatrixXd *m1, const C_Map_MatrixXd *m2, C_Map_MatrixXd *result)
  133. {
  134. c_to_eigen(result) = c_to_eigen(m1) * c_to_eigen(m2);
  135. }
  136. void Map_MatrixXd_add(const C_Map_MatrixXd *m1, const C_Map_MatrixXd *m2, C_Map_MatrixXd *result)
  137. {
  138. c_to_eigen(result) = c_to_eigen(m1) + c_to_eigen(m2);
  139. }