1/* matrix/gsl_matrix_char.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_CHAR_H__
21#define __GSL_MATRIX_CHAR_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_char.h>
29#include <gsl/gsl_blas_types.h>
30
31#undef __BEGIN_DECLS
32#undef __END_DECLS
33#ifdef __cplusplus
34# define __BEGIN_DECLS extern "C" {
35# define __END_DECLS }
36#else
37# define __BEGIN_DECLS /* empty */
38# define __END_DECLS /* empty */
39#endif
40
41__BEGIN_DECLS
42
43typedef struct
44{
45 size_t size1;
46 size_t size2;
47 size_t tda;
48 char * data;
49 gsl_block_char * block;
50 int owner;
51} gsl_matrix_char;
52
53typedef struct
54{
55 gsl_matrix_char matrix;
56} _gsl_matrix_char_view;
57
58typedef _gsl_matrix_char_view gsl_matrix_char_view;
59
60typedef struct
61{
62 gsl_matrix_char matrix;
63} _gsl_matrix_char_const_view;
64
65typedef const _gsl_matrix_char_const_view gsl_matrix_char_const_view;
66
67/* Allocation */
68
69gsl_matrix_char *
70gsl_matrix_char_alloc (const size_t n1, const size_t n2);
71
72gsl_matrix_char *
73gsl_matrix_char_calloc (const size_t n1, const size_t n2);
74
75gsl_matrix_char *
76gsl_matrix_char_alloc_from_block (gsl_block_char * b,
77 const size_t offset,
78 const size_t n1,
79 const size_t n2,
80 const size_t d2);
81
82gsl_matrix_char *
83gsl_matrix_char_alloc_from_matrix (gsl_matrix_char * m,
84 const size_t k1,
85 const size_t k2,
86 const size_t n1,
87 const size_t n2);
88
89gsl_vector_char *
90gsl_vector_char_alloc_row_from_matrix (gsl_matrix_char * m,
91 const size_t i);
92
93gsl_vector_char *
94gsl_vector_char_alloc_col_from_matrix (gsl_matrix_char * m,
95 const size_t j);
96
97void gsl_matrix_char_free (gsl_matrix_char * m);
98
99/* Views */
100
101_gsl_matrix_char_view
102gsl_matrix_char_submatrix (gsl_matrix_char * m,
103 const size_t i, const size_t j,
104 const size_t n1, const size_t n2);
105
106_gsl_vector_char_view
107gsl_matrix_char_row (gsl_matrix_char * m, const size_t i);
108
109_gsl_vector_char_view
110gsl_matrix_char_column (gsl_matrix_char * m, const size_t j);
111
112_gsl_vector_char_view
113gsl_matrix_char_diagonal (gsl_matrix_char * m);
114
115_gsl_vector_char_view
116gsl_matrix_char_subdiagonal (gsl_matrix_char * m, const size_t k);
117
118_gsl_vector_char_view
119gsl_matrix_char_superdiagonal (gsl_matrix_char * m, const size_t k);
120
121_gsl_vector_char_view
122gsl_matrix_char_subrow (gsl_matrix_char * m, const size_t i,
123 const size_t offset, const size_t n);
124
125_gsl_vector_char_view
126gsl_matrix_char_subcolumn (gsl_matrix_char * m, const size_t j,
127 const size_t offset, const size_t n);
128
129_gsl_matrix_char_view
130gsl_matrix_char_view_array (char * base,
131 const size_t n1,
132 const size_t n2);
133
134_gsl_matrix_char_view
135gsl_matrix_char_view_array_with_tda (char * base,
136 const size_t n1,
137 const size_t n2,
138 const size_t tda);
139
140
141_gsl_matrix_char_view
142gsl_matrix_char_view_vector (gsl_vector_char * v,
143 const size_t n1,
144 const size_t n2);
145
146_gsl_matrix_char_view
147gsl_matrix_char_view_vector_with_tda (gsl_vector_char * v,
148 const size_t n1,
149 const size_t n2,
150 const size_t tda);
151
152
153_gsl_matrix_char_const_view
154gsl_matrix_char_const_submatrix (const gsl_matrix_char * m,
155 const size_t i, const size_t j,
156 const size_t n1, const size_t n2);
157
158_gsl_vector_char_const_view
159gsl_matrix_char_const_row (const gsl_matrix_char * m,
160 const size_t i);
161
162_gsl_vector_char_const_view
163gsl_matrix_char_const_column (const gsl_matrix_char * m,
164 const size_t j);
165
166_gsl_vector_char_const_view
167gsl_matrix_char_const_diagonal (const gsl_matrix_char * m);
168
169_gsl_vector_char_const_view
170gsl_matrix_char_const_subdiagonal (const gsl_matrix_char * m,
171 const size_t k);
172
173_gsl_vector_char_const_view
174gsl_matrix_char_const_superdiagonal (const gsl_matrix_char * m,
175 const size_t k);
176
177_gsl_vector_char_const_view
178gsl_matrix_char_const_subrow (const gsl_matrix_char * m, const size_t i,
179 const size_t offset, const size_t n);
180
181_gsl_vector_char_const_view
182gsl_matrix_char_const_subcolumn (const gsl_matrix_char * m, const size_t j,
183 const size_t offset, const size_t n);
184
185_gsl_matrix_char_const_view
186gsl_matrix_char_const_view_array (const char * base,
187 const size_t n1,
188 const size_t n2);
189
190_gsl_matrix_char_const_view
191gsl_matrix_char_const_view_array_with_tda (const char * base,
192 const size_t n1,
193 const size_t n2,
194 const size_t tda);
195
196_gsl_matrix_char_const_view
197gsl_matrix_char_const_view_vector (const gsl_vector_char * v,
198 const size_t n1,
199 const size_t n2);
200
201_gsl_matrix_char_const_view
202gsl_matrix_char_const_view_vector_with_tda (const gsl_vector_char * v,
203 const size_t n1,
204 const size_t n2,
205 const size_t tda);
206
207/* Operations */
208
209void gsl_matrix_char_set_zero (gsl_matrix_char * m);
210void gsl_matrix_char_set_identity (gsl_matrix_char * m);
211void gsl_matrix_char_set_all (gsl_matrix_char * m, char x);
212
213int gsl_matrix_char_fread (FILE * stream, gsl_matrix_char * m) ;
214int gsl_matrix_char_fwrite (FILE * stream, const gsl_matrix_char * m) ;
215int gsl_matrix_char_fscanf (FILE * stream, gsl_matrix_char * m);
216int gsl_matrix_char_fprintf (FILE * stream, const gsl_matrix_char * m, const char * format);
217
218int gsl_matrix_char_memcpy(gsl_matrix_char * dest, const gsl_matrix_char * src);
219int gsl_matrix_char_swap(gsl_matrix_char * m1, gsl_matrix_char * m2);
220int gsl_matrix_char_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, gsl_matrix_char * dest, const gsl_matrix_char * src);
221
222int gsl_matrix_char_swap_rows(gsl_matrix_char * m, const size_t i, const size_t j);
223int gsl_matrix_char_swap_columns(gsl_matrix_char * m, const size_t i, const size_t j);
224int gsl_matrix_char_swap_rowcol(gsl_matrix_char * m, const size_t i, const size_t j);
225int gsl_matrix_char_transpose (gsl_matrix_char * m);
226int gsl_matrix_char_transpose_memcpy (gsl_matrix_char * dest, const gsl_matrix_char * src);
227int gsl_matrix_char_transpose_tricpy (CBLAS_UPLO_t Uplo_src, CBLAS_DIAG_t Diag, gsl_matrix_char * dest, const gsl_matrix_char * src);
228
229char gsl_matrix_char_max (const gsl_matrix_char * m);
230char gsl_matrix_char_min (const gsl_matrix_char * m);
231void gsl_matrix_char_minmax (const gsl_matrix_char * m, char * min_out, char * max_out);
232
233void gsl_matrix_char_max_index (const gsl_matrix_char * m, size_t * imax, size_t *jmax);
234void gsl_matrix_char_min_index (const gsl_matrix_char * m, size_t * imin, size_t *jmin);
235void gsl_matrix_char_minmax_index (const gsl_matrix_char * m, size_t * imin, size_t * jmin, size_t * imax, size_t * jmax);
236
237int gsl_matrix_char_equal (const gsl_matrix_char * a, const gsl_matrix_char * b);
238
239int gsl_matrix_char_isnull (const gsl_matrix_char * m);
240int gsl_matrix_char_ispos (const gsl_matrix_char * m);
241int gsl_matrix_char_isneg (const gsl_matrix_char * m);
242int gsl_matrix_char_isnonneg (const gsl_matrix_char * m);
243
244char gsl_matrix_char_norm1 (const gsl_matrix_char * m);
245
246int gsl_matrix_char_add (gsl_matrix_char * a, const gsl_matrix_char * b);
247int gsl_matrix_char_sub (gsl_matrix_char * a, const gsl_matrix_char * b);
248int gsl_matrix_char_mul_elements (gsl_matrix_char * a, const gsl_matrix_char * b);
249int gsl_matrix_char_div_elements (gsl_matrix_char * a, const gsl_matrix_char * b);
250int gsl_matrix_char_scale (gsl_matrix_char * a, const char x);
251int gsl_matrix_char_scale_rows (gsl_matrix_char * a, const gsl_vector_char * x);
252int gsl_matrix_char_scale_columns (gsl_matrix_char * a, const gsl_vector_char * x);
253int gsl_matrix_char_add_constant (gsl_matrix_char * a, const char x);
254int gsl_matrix_char_add_diagonal (gsl_matrix_char * a, const char x);
255
256/***********************************************************************/
257/* The functions below are obsolete */
258/***********************************************************************/
259int gsl_matrix_char_get_row(gsl_vector_char * v, const gsl_matrix_char * m, const size_t i);
260int gsl_matrix_char_get_col(gsl_vector_char * v, const gsl_matrix_char * m, const size_t j);
261int gsl_matrix_char_set_row(gsl_matrix_char * m, const size_t i, const gsl_vector_char * v);
262int gsl_matrix_char_set_col(gsl_matrix_char * m, const size_t j, const gsl_vector_char * v);
263/***********************************************************************/
264
265/* inline functions if you are using GCC */
266
267INLINE_DECL char gsl_matrix_char_get(const gsl_matrix_char * m, const size_t i, const size_t j);
268INLINE_DECL void gsl_matrix_char_set(gsl_matrix_char * m, const size_t i, const size_t j, const char x);
269INLINE_DECL char * gsl_matrix_char_ptr(gsl_matrix_char * m, const size_t i, const size_t j);
270INLINE_DECL const char * gsl_matrix_char_const_ptr(const gsl_matrix_char * m, const size_t i, const size_t j);
271
272#ifdef HAVE_INLINE
273INLINE_FUN
274char
275gsl_matrix_char_get(const gsl_matrix_char * m, const size_t i, const size_t j)
276{
277#if GSL_RANGE_CHECK
278 if (GSL_RANGE_COND(1))
279 {
280 if (i >= m->size1)
281 {
282 GSL_ERROR_VAL("first index out of range", GSL_EINVAL, 0) ;
283 }
284 else if (j >= m->size2)
285 {
286 GSL_ERROR_VAL("second index out of range", GSL_EINVAL, 0) ;
287 }
288 }
289#endif
290 return m->data[i * m->tda + j] ;
291}
292
293INLINE_FUN
294void
295gsl_matrix_char_set(gsl_matrix_char * m, const size_t i, const size_t j, const char x)
296{
297#if GSL_RANGE_CHECK
298 if (GSL_RANGE_COND(1))
299 {
300 if (i >= m->size1)
301 {
302 GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
303 }
304 else if (j >= m->size2)
305 {
306 GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
307 }
308 }
309#endif
310 m->data[i * m->tda + j] = x ;
311}
312
313INLINE_FUN
314char *
315gsl_matrix_char_ptr(gsl_matrix_char * m, const size_t i, const size_t j)
316{
317#if GSL_RANGE_CHECK
318 if (GSL_RANGE_COND(1))
319 {
320 if (i >= m->size1)
321 {
322 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
323 }
324 else if (j >= m->size2)
325 {
326 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
327 }
328 }
329#endif
330 return (char *) (m->data + (i * m->tda + j)) ;
331}
332
333INLINE_FUN
334const char *
335gsl_matrix_char_const_ptr(const gsl_matrix_char * m, const size_t i, const size_t j)
336{
337#if GSL_RANGE_CHECK
338 if (GSL_RANGE_COND(1))
339 {
340 if (i >= m->size1)
341 {
342 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
343 }
344 else if (j >= m->size2)
345 {
346 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
347 }
348 }
349#endif
350 return (const char *) (m->data + (i * m->tda + j)) ;
351}
352
353#endif
354
355__END_DECLS
356
357#endif /* __GSL_MATRIX_CHAR_H__ */
358