1 /* matrix/gsl_matrix_uint.h
3 * Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
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.
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.
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.
20 #ifndef __GSL_MATRIX_UINT_H__
21 #define __GSL_MATRIX_UINT_H__
24 #include <gsl/gsl_types.h>
25 #include <gsl/gsl_errno.h>
26 #include <gsl/gsl_check_range.h>
27 #include <gsl/gsl_vector_uint.h>
32 # define __BEGIN_DECLS extern "C" {
33 # define __END_DECLS }
35 # define __BEGIN_DECLS /* empty */
36 # define __END_DECLS /* empty */
47 gsl_block_uint * block;
53 gsl_matrix_uint matrix;
54 } _gsl_matrix_uint_view;
56 typedef _gsl_matrix_uint_view gsl_matrix_uint_view;
60 gsl_matrix_uint matrix;
61 } _gsl_matrix_uint_const_view;
63 typedef const _gsl_matrix_uint_const_view gsl_matrix_uint_const_view;
68 gsl_matrix_uint_alloc (const size_t n1, const size_t n2);
71 gsl_matrix_uint_calloc (const size_t n1, const size_t n2);
74 gsl_matrix_uint_alloc_from_block (gsl_block_uint * b,
81 gsl_matrix_uint_alloc_from_matrix (gsl_matrix_uint * m,
88 gsl_vector_uint_alloc_row_from_matrix (gsl_matrix_uint * m,
92 gsl_vector_uint_alloc_col_from_matrix (gsl_matrix_uint * m,
95 void gsl_matrix_uint_free (gsl_matrix_uint * m);
100 gsl_matrix_uint_submatrix (gsl_matrix_uint * m,
101 const size_t i, const size_t j,
102 const size_t n1, const size_t n2);
104 _gsl_vector_uint_view
105 gsl_matrix_uint_row (gsl_matrix_uint * m, const size_t i);
107 _gsl_vector_uint_view
108 gsl_matrix_uint_column (gsl_matrix_uint * m, const size_t j);
110 _gsl_vector_uint_view
111 gsl_matrix_uint_diagonal (gsl_matrix_uint * m);
113 _gsl_vector_uint_view
114 gsl_matrix_uint_subdiagonal (gsl_matrix_uint * m, const size_t k);
116 _gsl_vector_uint_view
117 gsl_matrix_uint_superdiagonal (gsl_matrix_uint * m, const size_t k);
119 _gsl_vector_uint_view
120 gsl_matrix_uint_subrow (gsl_matrix_uint * m, const size_t i,
121 const size_t offset, const size_t n);
123 _gsl_vector_uint_view
124 gsl_matrix_uint_subcolumn (gsl_matrix_uint * m, const size_t j,
125 const size_t offset, const size_t n);
127 _gsl_matrix_uint_view
128 gsl_matrix_uint_view_array (unsigned int * base,
132 _gsl_matrix_uint_view
133 gsl_matrix_uint_view_array_with_tda (unsigned int * base,
139 _gsl_matrix_uint_view
140 gsl_matrix_uint_view_vector (gsl_vector_uint * v,
144 _gsl_matrix_uint_view
145 gsl_matrix_uint_view_vector_with_tda (gsl_vector_uint * v,
151 _gsl_matrix_uint_const_view
152 gsl_matrix_uint_const_submatrix (const gsl_matrix_uint * m,
153 const size_t i, const size_t j,
154 const size_t n1, const size_t n2);
156 _gsl_vector_uint_const_view
157 gsl_matrix_uint_const_row (const gsl_matrix_uint * m,
160 _gsl_vector_uint_const_view
161 gsl_matrix_uint_const_column (const gsl_matrix_uint * m,
164 _gsl_vector_uint_const_view
165 gsl_matrix_uint_const_diagonal (const gsl_matrix_uint * m);
167 _gsl_vector_uint_const_view
168 gsl_matrix_uint_const_subdiagonal (const gsl_matrix_uint * m,
171 _gsl_vector_uint_const_view
172 gsl_matrix_uint_const_superdiagonal (const gsl_matrix_uint * m,
175 _gsl_vector_uint_const_view
176 gsl_matrix_uint_const_subrow (const gsl_matrix_uint * m, const size_t i,
177 const size_t offset, const size_t n);
179 _gsl_vector_uint_const_view
180 gsl_matrix_uint_const_subcolumn (const gsl_matrix_uint * m, const size_t j,
181 const size_t offset, const size_t n);
183 _gsl_matrix_uint_const_view
184 gsl_matrix_uint_const_view_array (const unsigned int * base,
188 _gsl_matrix_uint_const_view
189 gsl_matrix_uint_const_view_array_with_tda (const unsigned int * base,
194 _gsl_matrix_uint_const_view
195 gsl_matrix_uint_const_view_vector (const gsl_vector_uint * v,
199 _gsl_matrix_uint_const_view
200 gsl_matrix_uint_const_view_vector_with_tda (const gsl_vector_uint * v,
207 unsigned int gsl_matrix_uint_get(const gsl_matrix_uint * m, const size_t i, const size_t j);
208 void gsl_matrix_uint_set(gsl_matrix_uint * m, const size_t i, const size_t j, const unsigned int x);
210 unsigned int * gsl_matrix_uint_ptr(gsl_matrix_uint * m, const size_t i, const size_t j);
211 const unsigned int * gsl_matrix_uint_const_ptr(const gsl_matrix_uint * m, const size_t i, const size_t j);
213 void gsl_matrix_uint_set_zero (gsl_matrix_uint * m);
214 void gsl_matrix_uint_set_identity (gsl_matrix_uint * m);
215 void gsl_matrix_uint_set_all (gsl_matrix_uint * m, unsigned int x);
217 int gsl_matrix_uint_fread (FILE * stream, gsl_matrix_uint * m) ;
218 int gsl_matrix_uint_fwrite (FILE * stream, const gsl_matrix_uint * m) ;
219 int gsl_matrix_uint_fscanf (FILE * stream, gsl_matrix_uint * m);
220 int gsl_matrix_uint_fprintf (FILE * stream, const gsl_matrix_uint * m, const char * format);
222 int gsl_matrix_uint_memcpy(gsl_matrix_uint * dest, const gsl_matrix_uint * src);
223 int gsl_matrix_uint_swap(gsl_matrix_uint * m1, gsl_matrix_uint * m2);
225 int gsl_matrix_uint_swap_rows(gsl_matrix_uint * m, const size_t i, const size_t j);
226 int gsl_matrix_uint_swap_columns(gsl_matrix_uint * m, const size_t i, const size_t j);
227 int gsl_matrix_uint_swap_rowcol(gsl_matrix_uint * m, const size_t i, const size_t j);
228 int gsl_matrix_uint_transpose (gsl_matrix_uint * m);
229 int gsl_matrix_uint_transpose_memcpy (gsl_matrix_uint * dest, const gsl_matrix_uint * src);
231 unsigned int gsl_matrix_uint_max (const gsl_matrix_uint * m);
232 unsigned int gsl_matrix_uint_min (const gsl_matrix_uint * m);
233 void gsl_matrix_uint_minmax (const gsl_matrix_uint * m, unsigned int * min_out, unsigned int * max_out);
235 void gsl_matrix_uint_max_index (const gsl_matrix_uint * m, size_t * imax, size_t *jmax);
236 void gsl_matrix_uint_min_index (const gsl_matrix_uint * m, size_t * imin, size_t *jmin);
237 void gsl_matrix_uint_minmax_index (const gsl_matrix_uint * m, size_t * imin, size_t * jmin, size_t * imax, size_t * jmax);
239 int gsl_matrix_uint_isnull (const gsl_matrix_uint * m);
240 int gsl_matrix_uint_ispos (const gsl_matrix_uint * m);
241 int gsl_matrix_uint_isneg (const gsl_matrix_uint * m);
242 int gsl_matrix_uint_isnonneg (const gsl_matrix_uint * m);
244 int gsl_matrix_uint_add (gsl_matrix_uint * a, const gsl_matrix_uint * b);
245 int gsl_matrix_uint_sub (gsl_matrix_uint * a, const gsl_matrix_uint * b);
246 int gsl_matrix_uint_mul_elements (gsl_matrix_uint * a, const gsl_matrix_uint * b);
247 int gsl_matrix_uint_div_elements (gsl_matrix_uint * a, const gsl_matrix_uint * b);
248 int gsl_matrix_uint_scale (gsl_matrix_uint * a, const double x);
249 int gsl_matrix_uint_add_constant (gsl_matrix_uint * a, const double x);
250 int gsl_matrix_uint_add_diagonal (gsl_matrix_uint * a, const double x);
252 /***********************************************************************/
253 /* The functions below are obsolete */
254 /***********************************************************************/
255 int gsl_matrix_uint_get_row(gsl_vector_uint * v, const gsl_matrix_uint * m, const size_t i);
256 int gsl_matrix_uint_get_col(gsl_vector_uint * v, const gsl_matrix_uint * m, const size_t j);
257 int gsl_matrix_uint_set_row(gsl_matrix_uint * m, const size_t i, const gsl_vector_uint * v);
258 int gsl_matrix_uint_set_col(gsl_matrix_uint * m, const size_t j, const gsl_vector_uint * v);
260 /* inline functions if you are using GCC */
265 gsl_matrix_uint_get(const gsl_matrix_uint * m, const size_t i, const size_t j)
270 GSL_ERROR_VAL("first index out of range", GSL_EINVAL, 0) ;
272 else if (j >= m->size2)
274 GSL_ERROR_VAL("second index out of range", GSL_EINVAL, 0) ;
277 return m->data[i * m->tda + j] ;
282 gsl_matrix_uint_set(gsl_matrix_uint * m, const size_t i, const size_t j, const unsigned int x)
287 GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
289 else if (j >= m->size2)
291 GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
294 m->data[i * m->tda + j] = x ;
299 gsl_matrix_uint_ptr(gsl_matrix_uint * m, const size_t i, const size_t j)
304 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
306 else if (j >= m->size2)
308 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
311 return (unsigned int *) (m->data + (i * m->tda + j)) ;
316 gsl_matrix_uint_const_ptr(const gsl_matrix_uint * m, const size_t i, const size_t j)
321 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
323 else if (j >= m->size2)
325 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
328 return (const unsigned int *) (m->data + (i * m->tda + j)) ;
335 #endif /* __GSL_MATRIX_UINT_H__ */