gsl_matrix_double.h 10 KB

123456789101112131415161718192021222324252627282930313233343536373839404142434445464748495051525354555657585960616263646566676869707172737475767778798081828384858687888990919293949596979899100101102103104105106107108109110111112113114115116117118119120121122123124125126127128129130131132133134135136137138139140141142143144145146147148149150151152153154155156157158159160161162163164165166167168169170171172173174175176177178179180181182183184185186187188189190191192193194195196197198199200201202203204205206207208209210211212213214215216217218219220221222223224225226227228229230231232233234235236237238239240241242243244245246247248249250251252253254255256257258259260261262263264265266267268269270271272273274275276277278279280281282283284285286287288289290291292293294295296297298299300301302303304305306307308309310311312313314315316317318319320321322323324325326327328329330331332333334335336
  1. /* matrix/gsl_matrix_double.h
  2. *
  3. * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
  4. *
  5. * This program is free software; you can redistribute it and/or modify
  6. * it under the terms of the GNU General Public License as published by
  7. * the Free Software Foundation; either version 3 of the License, or (at
  8. * your option) any later version.
  9. *
  10. * This program is distributed in the hope that it will be useful, but
  11. * WITHOUT ANY WARRANTY; without even the implied warranty of
  12. * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
  13. * General Public License for more details.
  14. *
  15. * You should have received a copy of the GNU General Public License
  16. * along with this program; if not, write to the Free Software
  17. * Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
  18. */
  19. #ifndef __GSL_MATRIX_DOUBLE_H__
  20. #define __GSL_MATRIX_DOUBLE_H__
  21. #include <stdlib.h>
  22. #include "gsl_types.h"
  23. #include "gsl_errno.h"
  24. #include "gsl_check_range.h"
  25. #include "gsl_vector_double.h"
  26. #undef __BEGIN_DECLS
  27. #undef __END_DECLS
  28. #ifdef __cplusplus
  29. # define __BEGIN_DECLS extern "C" {
  30. # define __END_DECLS }
  31. #else
  32. # define __BEGIN_DECLS /* empty */
  33. # define __END_DECLS /* empty */
  34. #endif
  35. __BEGIN_DECLS
  36. typedef struct
  37. {
  38. size_t size1;
  39. size_t size2;
  40. size_t tda;
  41. double * data;
  42. gsl_block * block;
  43. int owner;
  44. } gsl_matrix;
  45. typedef struct
  46. {
  47. gsl_matrix matrix;
  48. } _gsl_matrix_view;
  49. typedef _gsl_matrix_view gsl_matrix_view;
  50. typedef struct
  51. {
  52. gsl_matrix matrix;
  53. } _gsl_matrix_const_view;
  54. typedef const _gsl_matrix_const_view gsl_matrix_const_view;
  55. /* Allocation */
  56. gsl_matrix *
  57. gsl_matrix_alloc (const size_t n1, const size_t n2);
  58. gsl_matrix *
  59. gsl_matrix_calloc (const size_t n1, const size_t n2);
  60. gsl_matrix *
  61. gsl_matrix_alloc_from_block (gsl_block * b,
  62. const size_t offset,
  63. const size_t n1,
  64. const size_t n2,
  65. const size_t d2);
  66. gsl_matrix *
  67. gsl_matrix_alloc_from_matrix (gsl_matrix * m,
  68. const size_t k1,
  69. const size_t k2,
  70. const size_t n1,
  71. const size_t n2);
  72. gsl_vector *
  73. gsl_vector_alloc_row_from_matrix (gsl_matrix * m,
  74. const size_t i);
  75. gsl_vector *
  76. gsl_vector_alloc_col_from_matrix (gsl_matrix * m,
  77. const size_t j);
  78. void gsl_matrix_free (gsl_matrix * m);
  79. /* Views */
  80. _gsl_matrix_view
  81. gsl_matrix_submatrix (gsl_matrix * m,
  82. const size_t i, const size_t j,
  83. const size_t n1, const size_t n2);
  84. _gsl_vector_view
  85. gsl_matrix_row (gsl_matrix * m, const size_t i);
  86. _gsl_vector_view
  87. gsl_matrix_column (gsl_matrix * m, const size_t j);
  88. _gsl_vector_view
  89. gsl_matrix_diagonal (gsl_matrix * m);
  90. _gsl_vector_view
  91. gsl_matrix_subdiagonal (gsl_matrix * m, const size_t k);
  92. _gsl_vector_view
  93. gsl_matrix_superdiagonal (gsl_matrix * m, const size_t k);
  94. _gsl_vector_view
  95. gsl_matrix_subrow (gsl_matrix * m, const size_t i,
  96. const size_t offset, const size_t n);
  97. _gsl_vector_view
  98. gsl_matrix_subcolumn (gsl_matrix * m, const size_t j,
  99. const size_t offset, const size_t n);
  100. _gsl_matrix_view
  101. gsl_matrix_view_array (double * base,
  102. const size_t n1,
  103. const size_t n2);
  104. _gsl_matrix_view
  105. gsl_matrix_view_array_with_tda (double * base,
  106. const size_t n1,
  107. const size_t n2,
  108. const size_t tda);
  109. _gsl_matrix_view
  110. gsl_matrix_view_vector (gsl_vector * v,
  111. const size_t n1,
  112. const size_t n2);
  113. _gsl_matrix_view
  114. gsl_matrix_view_vector_with_tda (gsl_vector * v,
  115. const size_t n1,
  116. const size_t n2,
  117. const size_t tda);
  118. _gsl_matrix_const_view
  119. gsl_matrix_const_submatrix (const gsl_matrix * m,
  120. const size_t i, const size_t j,
  121. const size_t n1, const size_t n2);
  122. _gsl_vector_const_view
  123. gsl_matrix_const_row (const gsl_matrix * m,
  124. const size_t i);
  125. _gsl_vector_const_view
  126. gsl_matrix_const_column (const gsl_matrix * m,
  127. const size_t j);
  128. _gsl_vector_const_view
  129. gsl_matrix_const_diagonal (const gsl_matrix * m);
  130. _gsl_vector_const_view
  131. gsl_matrix_const_subdiagonal (const gsl_matrix * m,
  132. const size_t k);
  133. _gsl_vector_const_view
  134. gsl_matrix_const_superdiagonal (const gsl_matrix * m,
  135. const size_t k);
  136. _gsl_vector_const_view
  137. gsl_matrix_const_subrow (const gsl_matrix * m, const size_t i,
  138. const size_t offset, const size_t n);
  139. _gsl_vector_const_view
  140. gsl_matrix_const_subcolumn (const gsl_matrix * m, const size_t j,
  141. const size_t offset, const size_t n);
  142. _gsl_matrix_const_view
  143. gsl_matrix_const_view_array (const double * base,
  144. const size_t n1,
  145. const size_t n2);
  146. _gsl_matrix_const_view
  147. gsl_matrix_const_view_array_with_tda (const double * base,
  148. const size_t n1,
  149. const size_t n2,
  150. const size_t tda);
  151. _gsl_matrix_const_view
  152. gsl_matrix_const_view_vector (const gsl_vector * v,
  153. const size_t n1,
  154. const size_t n2);
  155. _gsl_matrix_const_view
  156. gsl_matrix_const_view_vector_with_tda (const gsl_vector * v,
  157. const size_t n1,
  158. const size_t n2,
  159. const size_t tda);
  160. /* Operations */
  161. double gsl_matrix_get(const gsl_matrix * m, const size_t i, const size_t j);
  162. void gsl_matrix_set(gsl_matrix * m, const size_t i, const size_t j, const double x);
  163. double * gsl_matrix_ptr(gsl_matrix * m, const size_t i, const size_t j);
  164. const double * gsl_matrix_const_ptr(const gsl_matrix * m, const size_t i, const size_t j);
  165. void gsl_matrix_set_zero (gsl_matrix * m);
  166. void gsl_matrix_set_identity (gsl_matrix * m);
  167. void gsl_matrix_set_all (gsl_matrix * m, double x);
  168. int gsl_matrix_fread (FILE * stream, gsl_matrix * m) ;
  169. int gsl_matrix_fwrite (FILE * stream, const gsl_matrix * m) ;
  170. int gsl_matrix_fscanf (FILE * stream, gsl_matrix * m);
  171. int gsl_matrix_fprintf (FILE * stream, const gsl_matrix * m, const char * format);
  172. int gsl_matrix_memcpy(gsl_matrix * dest, const gsl_matrix * src);
  173. int gsl_matrix_swap(gsl_matrix * m1, gsl_matrix * m2);
  174. int gsl_matrix_swap_rows(gsl_matrix * m, const size_t i, const size_t j);
  175. int gsl_matrix_swap_columns(gsl_matrix * m, const size_t i, const size_t j);
  176. int gsl_matrix_swap_rowcol(gsl_matrix * m, const size_t i, const size_t j);
  177. int gsl_matrix_transpose (gsl_matrix * m);
  178. int gsl_matrix_transpose_memcpy (gsl_matrix * dest, const gsl_matrix * src);
  179. double gsl_matrix_max (const gsl_matrix * m);
  180. double gsl_matrix_min (const gsl_matrix * m);
  181. void gsl_matrix_minmax (const gsl_matrix * m, double * min_out, double * max_out);
  182. void gsl_matrix_max_index (const gsl_matrix * m, size_t * imax, size_t *jmax);
  183. void gsl_matrix_min_index (const gsl_matrix * m, size_t * imin, size_t *jmin);
  184. void gsl_matrix_minmax_index (const gsl_matrix * m, size_t * imin, size_t * jmin, size_t * imax, size_t * jmax);
  185. int gsl_matrix_isnull (const gsl_matrix * m);
  186. int gsl_matrix_ispos (const gsl_matrix * m);
  187. int gsl_matrix_isneg (const gsl_matrix * m);
  188. int gsl_matrix_isnonneg (const gsl_matrix * m);
  189. int gsl_matrix_add (gsl_matrix * a, const gsl_matrix * b);
  190. int gsl_matrix_sub (gsl_matrix * a, const gsl_matrix * b);
  191. int gsl_matrix_mul_elements (gsl_matrix * a, const gsl_matrix * b);
  192. int gsl_matrix_div_elements (gsl_matrix * a, const gsl_matrix * b);
  193. int gsl_matrix_scale (gsl_matrix * a, const double x);
  194. int gsl_matrix_add_constant (gsl_matrix * a, const double x);
  195. int gsl_matrix_add_diagonal (gsl_matrix * a, const double x);
  196. /***********************************************************************/
  197. /* The functions below are obsolete */
  198. /***********************************************************************/
  199. int gsl_matrix_get_row(gsl_vector * v, const gsl_matrix * m, const size_t i);
  200. int gsl_matrix_get_col(gsl_vector * v, const gsl_matrix * m, const size_t j);
  201. int gsl_matrix_set_row(gsl_matrix * m, const size_t i, const gsl_vector * v);
  202. int gsl_matrix_set_col(gsl_matrix * m, const size_t j, const gsl_vector * v);
  203. /* inline functions if you are using GCC */
  204. #ifdef HAVE_INLINE
  205. extern inline
  206. double
  207. gsl_matrix_get(const gsl_matrix * m, const size_t i, const size_t j)
  208. {
  209. #if GSL_RANGE_CHECK
  210. if (i >= m->size1)
  211. {
  212. GSL_ERROR_VAL("first index out of range", GSL_EINVAL, 0) ;
  213. }
  214. else if (j >= m->size2)
  215. {
  216. GSL_ERROR_VAL("second index out of range", GSL_EINVAL, 0) ;
  217. }
  218. #endif
  219. return m->data[i * m->tda + j] ;
  220. }
  221. extern inline
  222. void
  223. gsl_matrix_set(gsl_matrix * m, const size_t i, const size_t j, const double x)
  224. {
  225. #if GSL_RANGE_CHECK
  226. if (i >= m->size1)
  227. {
  228. GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
  229. }
  230. else if (j >= m->size2)
  231. {
  232. GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
  233. }
  234. #endif
  235. m->data[i * m->tda + j] = x ;
  236. }
  237. extern inline
  238. double *
  239. gsl_matrix_ptr(gsl_matrix * m, const size_t i, const size_t j)
  240. {
  241. #if GSL_RANGE_CHECK
  242. if (i >= m->size1)
  243. {
  244. GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
  245. }
  246. else if (j >= m->size2)
  247. {
  248. GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
  249. }
  250. #endif
  251. return (double *) (m->data + (i * m->tda + j)) ;
  252. }
  253. extern inline
  254. const double *
  255. gsl_matrix_const_ptr(const gsl_matrix * m, const size_t i, const size_t j)
  256. {
  257. #if GSL_RANGE_CHECK
  258. if (i >= m->size1)
  259. {
  260. GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
  261. }
  262. else if (j >= m->size2)
  263. {
  264. GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
  265. }
  266. #endif
  267. return (const double *) (m->data + (i * m->tda + j)) ;
  268. }
  269. #endif
  270. __END_DECLS
  271. #endif /* __GSL_MATRIX_DOUBLE_H__ */