1 /* matrix/gsl_matrix_int.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 
20 #ifndef __GSL_MATRIX_INT_H__
21 #define __GSL_MATRIX_INT_H__
22 
23 #include <stdlib.h>
24 #include <gsl/gsl_types.h>
25 #include <gsl/gsl_errno.h>
26 #include <gsl/gsl_inline.h>
27 #include <gsl/gsl_check_range.h>
28 #include <gsl/gsl_vector_int.h>
29 
30 #undef __BEGIN_DECLS
31 #undef __END_DECLS
32 #ifdef __cplusplus
33 # define __BEGIN_DECLS extern "C" {
34 # define __END_DECLS }
35 #else
36 # define __BEGIN_DECLS /* empty */
37 # define __END_DECLS /* empty */
38 #endif
39 
40 __BEGIN_DECLS
41 
42 typedef struct
43 {
44   size_t size1;
45   size_t size2;
46   size_t tda;
47   int * data;
48   gsl_block_int * block;
49   int owner;
50 } gsl_matrix_int;
51 
52 typedef struct
53 {
54   gsl_matrix_int matrix;
55 } _gsl_matrix_int_view;
56 
57 typedef _gsl_matrix_int_view gsl_matrix_int_view;
58 
59 typedef struct
60 {
61   gsl_matrix_int matrix;
62 } _gsl_matrix_int_const_view;
63 
64 typedef const _gsl_matrix_int_const_view gsl_matrix_int_const_view;
65 
66 /* Allocation */
67 
68 gsl_matrix_int *
69 gsl_matrix_int_alloc (const size_t n1, const size_t n2);
70 
71 gsl_matrix_int *
72 gsl_matrix_int_calloc (const size_t n1, const size_t n2);
73 
74 gsl_matrix_int *
75 gsl_matrix_int_alloc_from_block (gsl_block_int * b,
76                                    const size_t offset,
77                                    const size_t n1,
78                                    const size_t n2,
79                                    const size_t d2);
80 
81 gsl_matrix_int *
82 gsl_matrix_int_alloc_from_matrix (gsl_matrix_int * m,
83                                     const size_t k1,
84                                     const size_t k2,
85                                     const size_t n1,
86                                     const size_t n2);
87 
88 gsl_vector_int *
89 gsl_vector_int_alloc_row_from_matrix (gsl_matrix_int * m,
90                                         const size_t i);
91 
92 gsl_vector_int *
93 gsl_vector_int_alloc_col_from_matrix (gsl_matrix_int * m,
94                                         const size_t j);
95 
96 void gsl_matrix_int_free (gsl_matrix_int * m);
97 
98 /* Views */
99 
100 _gsl_matrix_int_view
101 gsl_matrix_int_submatrix (gsl_matrix_int * m,
102                             const size_t i, const size_t j,
103                             const size_t n1, const size_t n2);
104 
105 _gsl_vector_int_view
106 gsl_matrix_int_row (gsl_matrix_int * m, const size_t i);
107 
108 _gsl_vector_int_view
109 gsl_matrix_int_column (gsl_matrix_int * m, const size_t j);
110 
111 _gsl_vector_int_view
112 gsl_matrix_int_diagonal (gsl_matrix_int * m);
113 
114 _gsl_vector_int_view
115 gsl_matrix_int_subdiagonal (gsl_matrix_int * m, const size_t k);
116 
117 _gsl_vector_int_view
118 gsl_matrix_int_superdiagonal (gsl_matrix_int * m, const size_t k);
119 
120 _gsl_vector_int_view
121 gsl_matrix_int_subrow (gsl_matrix_int * m, const size_t i,
122                          const size_t offset, const size_t n);
123 
124 _gsl_vector_int_view
125 gsl_matrix_int_subcolumn (gsl_matrix_int * m, const size_t j,
126                             const size_t offset, const size_t n);
127 
128 _gsl_matrix_int_view
129 gsl_matrix_int_view_array (int * base,
130                              const size_t n1,
131                              const size_t n2);
132 
133 _gsl_matrix_int_view
134 gsl_matrix_int_view_array_with_tda (int * base,
135                                       const size_t n1,
136                                       const size_t n2,
137                                       const size_t tda);
138 
139 
140 _gsl_matrix_int_view
141 gsl_matrix_int_view_vector (gsl_vector_int * v,
142                               const size_t n1,
143                               const size_t n2);
144 
145 _gsl_matrix_int_view
146 gsl_matrix_int_view_vector_with_tda (gsl_vector_int * v,
147                                        const size_t n1,
148                                        const size_t n2,
149                                        const size_t tda);
150 
151 
152 _gsl_matrix_int_const_view
153 gsl_matrix_int_const_submatrix (const gsl_matrix_int * m,
154                                   const size_t i, const size_t j,
155                                   const size_t n1, const size_t n2);
156 
157 _gsl_vector_int_const_view
158 gsl_matrix_int_const_row (const gsl_matrix_int * m,
159                             const size_t i);
160 
161 _gsl_vector_int_const_view
162 gsl_matrix_int_const_column (const gsl_matrix_int * m,
163                                const size_t j);
164 
165 _gsl_vector_int_const_view
166 gsl_matrix_int_const_diagonal (const gsl_matrix_int * m);
167 
168 _gsl_vector_int_const_view
169 gsl_matrix_int_const_subdiagonal (const gsl_matrix_int * m,
170                                     const size_t k);
171 
172 _gsl_vector_int_const_view
173 gsl_matrix_int_const_superdiagonal (const gsl_matrix_int * m,
174                                       const size_t k);
175 
176 _gsl_vector_int_const_view
177 gsl_matrix_int_const_subrow (const gsl_matrix_int * m, const size_t i,
178                                const size_t offset, const size_t n);
179 
180 _gsl_vector_int_const_view
181 gsl_matrix_int_const_subcolumn (const gsl_matrix_int * m, const size_t j,
182                                   const size_t offset, const size_t n);
183 
184 _gsl_matrix_int_const_view
185 gsl_matrix_int_const_view_array (const int * base,
186                                    const size_t n1,
187                                    const size_t n2);
188 
189 _gsl_matrix_int_const_view
190 gsl_matrix_int_const_view_array_with_tda (const int * base,
191                                             const size_t n1,
192                                             const size_t n2,
193                                             const size_t tda);
194 
195 _gsl_matrix_int_const_view
196 gsl_matrix_int_const_view_vector (const gsl_vector_int * v,
197                                     const size_t n1,
198                                     const size_t n2);
199 
200 _gsl_matrix_int_const_view
201 gsl_matrix_int_const_view_vector_with_tda (const gsl_vector_int * v,
202                                              const size_t n1,
203                                              const size_t n2,
204                                              const size_t tda);
205 
206 /* Operations */
207 
208 void gsl_matrix_int_set_zero (gsl_matrix_int * m);
209 void gsl_matrix_int_set_identity (gsl_matrix_int * m);
210 void gsl_matrix_int_set_all (gsl_matrix_int * m, int x);
211 
212 int gsl_matrix_int_fread (FILE * stream, gsl_matrix_int * m) ;
213 int gsl_matrix_int_fwrite (FILE * stream, const gsl_matrix_int * m) ;
214 int gsl_matrix_int_fscanf (FILE * stream, gsl_matrix_int * m);
215 int gsl_matrix_int_fprintf (FILE * stream, const gsl_matrix_int * m, const char * format);
216 
217 int gsl_matrix_int_memcpy(gsl_matrix_int * dest, const gsl_matrix_int * src);
218 int gsl_matrix_int_swap(gsl_matrix_int * m1, gsl_matrix_int * m2);
219 
220 int gsl_matrix_int_swap_rows(gsl_matrix_int * m, const size_t i, const size_t j);
221 int gsl_matrix_int_swap_columns(gsl_matrix_int * m, const size_t i, const size_t j);
222 int gsl_matrix_int_swap_rowcol(gsl_matrix_int * m, const size_t i, const size_t j);
223 int gsl_matrix_int_transpose (gsl_matrix_int * m);
224 int gsl_matrix_int_transpose_memcpy (gsl_matrix_int * dest, const gsl_matrix_int * src);
225 
226 int gsl_matrix_int_max (const gsl_matrix_int * m);
227 int gsl_matrix_int_min (const gsl_matrix_int * m);
228 void gsl_matrix_int_minmax (const gsl_matrix_int * m, int * min_out, int * max_out);
229 
230 void gsl_matrix_int_max_index (const gsl_matrix_int * m, size_t * imax, size_t *jmax);
231 void gsl_matrix_int_min_index (const gsl_matrix_int * m, size_t * imin, size_t *jmin);
232 void gsl_matrix_int_minmax_index (const gsl_matrix_int * m, size_t * imin, size_t * jmin, size_t * imax, size_t * jmax);
233 
234 int gsl_matrix_int_equal (const gsl_matrix_int * a, const gsl_matrix_int * b);
235 
236 int gsl_matrix_int_isnull (const gsl_matrix_int * m);
237 int gsl_matrix_int_ispos (const gsl_matrix_int * m);
238 int gsl_matrix_int_isneg (const gsl_matrix_int * m);
239 int gsl_matrix_int_isnonneg (const gsl_matrix_int * m);
240 
241 int gsl_matrix_int_add (gsl_matrix_int * a, const gsl_matrix_int * b);
242 int gsl_matrix_int_sub (gsl_matrix_int * a, const gsl_matrix_int * b);
243 int gsl_matrix_int_mul_elements (gsl_matrix_int * a, const gsl_matrix_int * b);
244 int gsl_matrix_int_div_elements (gsl_matrix_int * a, const gsl_matrix_int * b);
245 int gsl_matrix_int_scale (gsl_matrix_int * a, const double x);
246 int gsl_matrix_int_add_constant (gsl_matrix_int * a, const double x);
247 int gsl_matrix_int_add_diagonal (gsl_matrix_int * a, const double x);
248 
249 /***********************************************************************/
250 /* The functions below are obsolete                                    */
251 /***********************************************************************/
252 int gsl_matrix_int_get_row(gsl_vector_int * v, const gsl_matrix_int * m, const size_t i);
253 int gsl_matrix_int_get_col(gsl_vector_int * v, const gsl_matrix_int * m, const size_t j);
254 int gsl_matrix_int_set_row(gsl_matrix_int * m, const size_t i, const gsl_vector_int * v);
255 int gsl_matrix_int_set_col(gsl_matrix_int * m, const size_t j, const gsl_vector_int * v);
256 /***********************************************************************/
257 
258 /* inline functions if you are using GCC */
259 
260 INLINE_DECL int   gsl_matrix_int_get(const gsl_matrix_int * m, const size_t i, const size_t j);
261 INLINE_DECL void    gsl_matrix_int_set(gsl_matrix_int * m, const size_t i, const size_t j, const int x);
262 INLINE_DECL int * gsl_matrix_int_ptr(gsl_matrix_int * m, const size_t i, const size_t j);
263 INLINE_DECL const int * gsl_matrix_int_const_ptr(const gsl_matrix_int * m, const size_t i, const size_t j);
264 
265 #ifdef HAVE_INLINE
266 INLINE_FUN
267 int
gsl_matrix_int_get(const gsl_matrix_int * m,const size_t i,const size_t j)268 gsl_matrix_int_get(const gsl_matrix_int * m, const size_t i, const size_t j)
269 {
270 #if GSL_RANGE_CHECK
271   if (GSL_RANGE_COND(1))
272     {
273       if (i >= m->size1)
274         {
275           GSL_ERROR_VAL("first index out of range", GSL_EINVAL, 0) ;
276         }
277       else if (j >= m->size2)
278         {
279           GSL_ERROR_VAL("second index out of range", GSL_EINVAL, 0) ;
280         }
281     }
282 #endif
283   return m->data[i * m->tda + j] ;
284 }
285 
286 INLINE_FUN
287 void
gsl_matrix_int_set(gsl_matrix_int * m,const size_t i,const size_t j,const int x)288 gsl_matrix_int_set(gsl_matrix_int * m, const size_t i, const size_t j, const int x)
289 {
290 #if GSL_RANGE_CHECK
291   if (GSL_RANGE_COND(1))
292     {
293       if (i >= m->size1)
294         {
295           GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
296         }
297       else if (j >= m->size2)
298         {
299           GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
300         }
301     }
302 #endif
303   m->data[i * m->tda + j] = x ;
304 }
305 
306 INLINE_FUN
307 int *
gsl_matrix_int_ptr(gsl_matrix_int * m,const size_t i,const size_t j)308 gsl_matrix_int_ptr(gsl_matrix_int * m, const size_t i, const size_t j)
309 {
310 #if GSL_RANGE_CHECK
311   if (GSL_RANGE_COND(1))
312     {
313       if (i >= m->size1)
314         {
315           GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
316         }
317       else if (j >= m->size2)
318         {
319           GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
320         }
321     }
322 #endif
323   return (int *) (m->data + (i * m->tda + j)) ;
324 }
325 
326 INLINE_FUN
327 const int *
gsl_matrix_int_const_ptr(const gsl_matrix_int * m,const size_t i,const size_t j)328 gsl_matrix_int_const_ptr(const gsl_matrix_int * m, const size_t i, const size_t j)
329 {
330 #if GSL_RANGE_CHECK
331   if (GSL_RANGE_COND(1))
332     {
333       if (i >= m->size1)
334         {
335           GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
336         }
337       else if (j >= m->size2)
338         {
339           GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
340         }
341     }
342 #endif
343   return (const int *) (m->data + (i * m->tda + j)) ;
344 }
345 
346 #endif
347 
348 __END_DECLS
349 
350 #endif /* __GSL_MATRIX_INT_H__ */
351