1 /* matrix/gsl_matrix_complex_double.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_COMPLEX_DOUBLE_H__
21 #define __GSL_MATRIX_COMPLEX_DOUBLE_H__
24 #include <gsl/gsl_types.h>
25 #include <gsl/gsl_errno.h>
26 #include <gsl/gsl_complex.h>
27 #include <gsl/gsl_check_range.h>
28 #include <gsl/gsl_vector_complex_double.h>
29 #include <gsl/gsl_blas_types.h>
34 # define __BEGIN_DECLS extern "C" {
35 # define __END_DECLS }
37 # define __BEGIN_DECLS /* empty */
38 # define __END_DECLS /* empty */
49 gsl_block_complex * block;
51 } gsl_matrix_complex ;
55 gsl_matrix_complex matrix;
56 } _gsl_matrix_complex_view;
58 typedef _gsl_matrix_complex_view gsl_matrix_complex_view;
62 gsl_matrix_complex matrix;
63 } _gsl_matrix_complex_const_view;
65 typedef const _gsl_matrix_complex_const_view gsl_matrix_complex_const_view;
71 gsl_matrix_complex_alloc (const size_t n1, const size_t n2);
74 gsl_matrix_complex_calloc (const size_t n1, const size_t n2);
77 gsl_matrix_complex_alloc_from_block (gsl_block_complex * b,
79 const size_t n1, const size_t n2, const size_t d2);
82 gsl_matrix_complex_alloc_from_matrix (gsl_matrix_complex * b,
83 const size_t k1, const size_t k2,
84 const size_t n1, const size_t n2);
87 gsl_vector_complex_alloc_row_from_matrix (gsl_matrix_complex * m,
91 gsl_vector_complex_alloc_col_from_matrix (gsl_matrix_complex * m,
94 void gsl_matrix_complex_free (gsl_matrix_complex * m);
98 _gsl_matrix_complex_view
99 gsl_matrix_complex_submatrix (gsl_matrix_complex * m,
100 const size_t i, const size_t j,
101 const size_t n1, const size_t n2);
103 _gsl_vector_complex_view
104 gsl_matrix_complex_row (gsl_matrix_complex * m, const size_t i);
106 _gsl_vector_complex_view
107 gsl_matrix_complex_column (gsl_matrix_complex * m, const size_t j);
109 _gsl_vector_complex_view
110 gsl_matrix_complex_diagonal (gsl_matrix_complex * m);
112 _gsl_vector_complex_view
113 gsl_matrix_complex_subdiagonal (gsl_matrix_complex * m, const size_t k);
115 _gsl_vector_complex_view
116 gsl_matrix_complex_superdiagonal (gsl_matrix_complex * m, const size_t k);
118 _gsl_vector_complex_view
119 gsl_matrix_complex_subrow (gsl_matrix_complex * m,
120 const size_t i, const size_t offset,
123 _gsl_vector_complex_view
124 gsl_matrix_complex_subcolumn (gsl_matrix_complex * m,
125 const size_t j, const size_t offset,
128 _gsl_matrix_complex_view
129 gsl_matrix_complex_view_array (double * base,
133 _gsl_matrix_complex_view
134 gsl_matrix_complex_view_array_with_tda (double * base,
139 _gsl_matrix_complex_view
140 gsl_matrix_complex_view_vector (gsl_vector_complex * v,
144 _gsl_matrix_complex_view
145 gsl_matrix_complex_view_vector_with_tda (gsl_vector_complex * v,
151 _gsl_matrix_complex_const_view
152 gsl_matrix_complex_const_submatrix (const gsl_matrix_complex * m,
153 const size_t i, const size_t j,
154 const size_t n1, const size_t n2);
156 _gsl_vector_complex_const_view
157 gsl_matrix_complex_const_row (const gsl_matrix_complex * m,
160 _gsl_vector_complex_const_view
161 gsl_matrix_complex_const_column (const gsl_matrix_complex * m,
164 _gsl_vector_complex_const_view
165 gsl_matrix_complex_const_diagonal (const gsl_matrix_complex * m);
167 _gsl_vector_complex_const_view
168 gsl_matrix_complex_const_subdiagonal (const gsl_matrix_complex * m,
171 _gsl_vector_complex_const_view
172 gsl_matrix_complex_const_superdiagonal (const gsl_matrix_complex * m,
175 _gsl_vector_complex_const_view
176 gsl_matrix_complex_const_subrow (const gsl_matrix_complex * m,
177 const size_t i, const size_t offset,
180 _gsl_vector_complex_const_view
181 gsl_matrix_complex_const_subcolumn (const gsl_matrix_complex * m,
182 const size_t j, const size_t offset,
185 _gsl_matrix_complex_const_view
186 gsl_matrix_complex_const_view_array (const double * base,
190 _gsl_matrix_complex_const_view
191 gsl_matrix_complex_const_view_array_with_tda (const double * base,
196 _gsl_matrix_complex_const_view
197 gsl_matrix_complex_const_view_vector (const gsl_vector_complex * v,
201 _gsl_matrix_complex_const_view
202 gsl_matrix_complex_const_view_vector_with_tda (const gsl_vector_complex * v,
209 void gsl_matrix_complex_set_zero (gsl_matrix_complex * m);
210 void gsl_matrix_complex_set_identity (gsl_matrix_complex * m);
211 void gsl_matrix_complex_set_all (gsl_matrix_complex * m, gsl_complex x);
213 int gsl_matrix_complex_fread (FILE * stream, gsl_matrix_complex * m) ;
214 int gsl_matrix_complex_fwrite (FILE * stream, const gsl_matrix_complex * m) ;
215 int gsl_matrix_complex_fscanf (FILE * stream, gsl_matrix_complex * m);
216 int gsl_matrix_complex_fprintf (FILE * stream, const gsl_matrix_complex * m, const char * format);
218 int gsl_matrix_complex_memcpy(gsl_matrix_complex * dest, const gsl_matrix_complex * src);
219 int gsl_matrix_complex_swap(gsl_matrix_complex * m1, gsl_matrix_complex * m2);
220 int gsl_matrix_complex_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, gsl_matrix_complex * dest, const gsl_matrix_complex * src);
222 int gsl_matrix_complex_swap_rows(gsl_matrix_complex * m, const size_t i, const size_t j);
223 int gsl_matrix_complex_swap_columns(gsl_matrix_complex * m, const size_t i, const size_t j);
224 int gsl_matrix_complex_swap_rowcol(gsl_matrix_complex * m, const size_t i, const size_t j);
226 int gsl_matrix_complex_transpose (gsl_matrix_complex * m);
227 int gsl_matrix_complex_transpose_memcpy (gsl_matrix_complex * dest, const gsl_matrix_complex * src);
228 int gsl_matrix_complex_transpose_tricpy(CBLAS_UPLO_t Uplo_src, CBLAS_DIAG_t Diag, gsl_matrix_complex * dest, const gsl_matrix_complex * src);
230 int gsl_matrix_complex_conjtrans_memcpy (gsl_matrix_complex * dest, const gsl_matrix_complex * src);
232 int gsl_matrix_complex_equal (const gsl_matrix_complex * a, const gsl_matrix_complex * b);
234 int gsl_matrix_complex_isnull (const gsl_matrix_complex * m);
235 int gsl_matrix_complex_ispos (const gsl_matrix_complex * m);
236 int gsl_matrix_complex_isneg (const gsl_matrix_complex * m);
237 int gsl_matrix_complex_isnonneg (const gsl_matrix_complex * m);
239 int gsl_matrix_complex_add (gsl_matrix_complex * a, const gsl_matrix_complex * b);
240 int gsl_matrix_complex_sub (gsl_matrix_complex * a, const gsl_matrix_complex * b);
241 int gsl_matrix_complex_mul_elements (gsl_matrix_complex * a, const gsl_matrix_complex * b);
242 int gsl_matrix_complex_div_elements (gsl_matrix_complex * a, const gsl_matrix_complex * b);
243 int gsl_matrix_complex_scale (gsl_matrix_complex * a, const gsl_complex x);
244 int gsl_matrix_complex_scale_rows (gsl_matrix_complex * a, const gsl_vector_complex * x);
245 int gsl_matrix_complex_scale_columns (gsl_matrix_complex * a, const gsl_vector_complex * x);
246 int gsl_matrix_complex_add_constant (gsl_matrix_complex * a, const gsl_complex x);
247 int gsl_matrix_complex_add_diagonal (gsl_matrix_complex * a, const gsl_complex x);
249 /***********************************************************************/
250 /* The functions below are obsolete */
251 /***********************************************************************/
252 int gsl_matrix_complex_get_row(gsl_vector_complex * v, const gsl_matrix_complex * m, const size_t i);
253 int gsl_matrix_complex_get_col(gsl_vector_complex * v, const gsl_matrix_complex * m, const size_t j);
254 int gsl_matrix_complex_set_row(gsl_matrix_complex * m, const size_t i, const gsl_vector_complex * v);
255 int gsl_matrix_complex_set_col(gsl_matrix_complex * m, const size_t j, const gsl_vector_complex * v);
256 /***********************************************************************/
258 /* inline functions if you are using GCC */
260 INLINE_DECL gsl_complex gsl_matrix_complex_get(const gsl_matrix_complex * m, const size_t i, const size_t j);
261 INLINE_DECL void gsl_matrix_complex_set(gsl_matrix_complex * m, const size_t i, const size_t j, const gsl_complex x);
263 INLINE_DECL gsl_complex * gsl_matrix_complex_ptr(gsl_matrix_complex * m, const size_t i, const size_t j);
264 INLINE_DECL const gsl_complex * gsl_matrix_complex_const_ptr(const gsl_matrix_complex * m, const size_t i, const size_t j);
270 gsl_matrix_complex_get(const gsl_matrix_complex * m,
271 const size_t i, const size_t j)
274 if (GSL_RANGE_COND(1))
276 gsl_complex zero = {{0,0}};
280 GSL_ERROR_VAL("first index out of range", GSL_EINVAL, zero) ;
282 else if (j >= m->size2)
284 GSL_ERROR_VAL("second index out of range", GSL_EINVAL, zero) ;
288 return *(gsl_complex *)(m->data + 2*(i * m->tda + j)) ;
293 gsl_matrix_complex_set(gsl_matrix_complex * m,
294 const size_t i, const size_t j, const gsl_complex x)
297 if (GSL_RANGE_COND(1))
301 GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
303 else if (j >= m->size2)
305 GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
309 *(gsl_complex *)(m->data + 2*(i * m->tda + j)) = x ;
314 gsl_matrix_complex_ptr(gsl_matrix_complex * m,
315 const size_t i, const size_t j)
318 if (GSL_RANGE_COND(1))
322 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
324 else if (j >= m->size2)
326 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
330 return (gsl_complex *)(m->data + 2*(i * m->tda + j)) ;
335 gsl_matrix_complex_const_ptr(const gsl_matrix_complex * m,
336 const size_t i, const size_t j)
339 if (GSL_RANGE_COND(1))
343 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
345 else if (j >= m->size2)
347 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
351 return (const gsl_complex *)(m->data + 2*(i * m->tda + j)) ;
354 #endif /* HAVE_INLINE */
358 #endif /* __GSL_MATRIX_COMPLEX_DOUBLE_H__ */