File indexing completed on 2025-02-21 10:03:49
0001
0002
0003
0004
0005
0006
0007
0008
0009
0010
0011
0012
0013
0014
0015
0016
0017
0018
0019
0020 #ifndef __GSL_MATRIX_COMPLEX_DOUBLE_H__
0021 #define __GSL_MATRIX_COMPLEX_DOUBLE_H__
0022
0023 #include <stdlib.h>
0024 #include <gsl/gsl_types.h>
0025 #include <gsl/gsl_errno.h>
0026 #include <gsl/gsl_complex.h>
0027 #include <gsl/gsl_check_range.h>
0028 #include <gsl/gsl_vector_complex_double.h>
0029 #include <gsl/gsl_blas_types.h>
0030
0031 #undef __BEGIN_DECLS
0032 #undef __END_DECLS
0033 #ifdef __cplusplus
0034 # define __BEGIN_DECLS extern "C" {
0035 # define __END_DECLS }
0036 #else
0037 # define __BEGIN_DECLS
0038 # define __END_DECLS
0039 #endif
0040
0041 __BEGIN_DECLS
0042
0043 typedef struct
0044 {
0045 size_t size1;
0046 size_t size2;
0047 size_t tda;
0048 double * data;
0049 gsl_block_complex * block;
0050 int owner;
0051 } gsl_matrix_complex ;
0052
0053 typedef struct
0054 {
0055 gsl_matrix_complex matrix;
0056 } _gsl_matrix_complex_view;
0057
0058 typedef _gsl_matrix_complex_view gsl_matrix_complex_view;
0059
0060 typedef struct
0061 {
0062 gsl_matrix_complex matrix;
0063 } _gsl_matrix_complex_const_view;
0064
0065 typedef const _gsl_matrix_complex_const_view gsl_matrix_complex_const_view;
0066
0067
0068
0069
0070 gsl_matrix_complex *
0071 gsl_matrix_complex_alloc (const size_t n1, const size_t n2);
0072
0073 gsl_matrix_complex *
0074 gsl_matrix_complex_calloc (const size_t n1, const size_t n2);
0075
0076 gsl_matrix_complex *
0077 gsl_matrix_complex_alloc_from_block (gsl_block_complex * b,
0078 const size_t offset,
0079 const size_t n1, const size_t n2, const size_t d2);
0080
0081 gsl_matrix_complex *
0082 gsl_matrix_complex_alloc_from_matrix (gsl_matrix_complex * b,
0083 const size_t k1, const size_t k2,
0084 const size_t n1, const size_t n2);
0085
0086 gsl_vector_complex *
0087 gsl_vector_complex_alloc_row_from_matrix (gsl_matrix_complex * m,
0088 const size_t i);
0089
0090 gsl_vector_complex *
0091 gsl_vector_complex_alloc_col_from_matrix (gsl_matrix_complex * m,
0092 const size_t j);
0093
0094 void gsl_matrix_complex_free (gsl_matrix_complex * m);
0095
0096
0097
0098 _gsl_matrix_complex_view
0099 gsl_matrix_complex_submatrix (gsl_matrix_complex * m,
0100 const size_t i, const size_t j,
0101 const size_t n1, const size_t n2);
0102
0103 _gsl_vector_complex_view
0104 gsl_matrix_complex_row (gsl_matrix_complex * m, const size_t i);
0105
0106 _gsl_vector_complex_view
0107 gsl_matrix_complex_column (gsl_matrix_complex * m, const size_t j);
0108
0109 _gsl_vector_complex_view
0110 gsl_matrix_complex_diagonal (gsl_matrix_complex * m);
0111
0112 _gsl_vector_complex_view
0113 gsl_matrix_complex_subdiagonal (gsl_matrix_complex * m, const size_t k);
0114
0115 _gsl_vector_complex_view
0116 gsl_matrix_complex_superdiagonal (gsl_matrix_complex * m, const size_t k);
0117
0118 _gsl_vector_complex_view
0119 gsl_matrix_complex_subrow (gsl_matrix_complex * m,
0120 const size_t i, const size_t offset,
0121 const size_t n);
0122
0123 _gsl_vector_complex_view
0124 gsl_matrix_complex_subcolumn (gsl_matrix_complex * m,
0125 const size_t j, const size_t offset,
0126 const size_t n);
0127
0128 _gsl_matrix_complex_view
0129 gsl_matrix_complex_view_array (double * base,
0130 const size_t n1,
0131 const size_t n2);
0132
0133 _gsl_matrix_complex_view
0134 gsl_matrix_complex_view_array_with_tda (double * base,
0135 const size_t n1,
0136 const size_t n2,
0137 const size_t tda);
0138
0139 _gsl_matrix_complex_view
0140 gsl_matrix_complex_view_vector (gsl_vector_complex * v,
0141 const size_t n1,
0142 const size_t n2);
0143
0144 _gsl_matrix_complex_view
0145 gsl_matrix_complex_view_vector_with_tda (gsl_vector_complex * v,
0146 const size_t n1,
0147 const size_t n2,
0148 const size_t tda);
0149
0150
0151 _gsl_matrix_complex_const_view
0152 gsl_matrix_complex_const_submatrix (const gsl_matrix_complex * m,
0153 const size_t i, const size_t j,
0154 const size_t n1, const size_t n2);
0155
0156 _gsl_vector_complex_const_view
0157 gsl_matrix_complex_const_row (const gsl_matrix_complex * m,
0158 const size_t i);
0159
0160 _gsl_vector_complex_const_view
0161 gsl_matrix_complex_const_column (const gsl_matrix_complex * m,
0162 const size_t j);
0163
0164 _gsl_vector_complex_const_view
0165 gsl_matrix_complex_const_diagonal (const gsl_matrix_complex * m);
0166
0167 _gsl_vector_complex_const_view
0168 gsl_matrix_complex_const_subdiagonal (const gsl_matrix_complex * m,
0169 const size_t k);
0170
0171 _gsl_vector_complex_const_view
0172 gsl_matrix_complex_const_superdiagonal (const gsl_matrix_complex * m,
0173 const size_t k);
0174
0175 _gsl_vector_complex_const_view
0176 gsl_matrix_complex_const_subrow (const gsl_matrix_complex * m,
0177 const size_t i, const size_t offset,
0178 const size_t n);
0179
0180 _gsl_vector_complex_const_view
0181 gsl_matrix_complex_const_subcolumn (const gsl_matrix_complex * m,
0182 const size_t j, const size_t offset,
0183 const size_t n);
0184
0185 _gsl_matrix_complex_const_view
0186 gsl_matrix_complex_const_view_array (const double * base,
0187 const size_t n1,
0188 const size_t n2);
0189
0190 _gsl_matrix_complex_const_view
0191 gsl_matrix_complex_const_view_array_with_tda (const double * base,
0192 const size_t n1,
0193 const size_t n2,
0194 const size_t tda);
0195
0196 _gsl_matrix_complex_const_view
0197 gsl_matrix_complex_const_view_vector (const gsl_vector_complex * v,
0198 const size_t n1,
0199 const size_t n2);
0200
0201 _gsl_matrix_complex_const_view
0202 gsl_matrix_complex_const_view_vector_with_tda (const gsl_vector_complex * v,
0203 const size_t n1,
0204 const size_t n2,
0205 const size_t tda);
0206
0207
0208
0209 void gsl_matrix_complex_set_zero (gsl_matrix_complex * m);
0210 void gsl_matrix_complex_set_identity (gsl_matrix_complex * m);
0211 void gsl_matrix_complex_set_all (gsl_matrix_complex * m, gsl_complex x);
0212
0213 int gsl_matrix_complex_fread (FILE * stream, gsl_matrix_complex * m) ;
0214 int gsl_matrix_complex_fwrite (FILE * stream, const gsl_matrix_complex * m) ;
0215 int gsl_matrix_complex_fscanf (FILE * stream, gsl_matrix_complex * m);
0216 int gsl_matrix_complex_fprintf (FILE * stream, const gsl_matrix_complex * m, const char * format);
0217
0218 int gsl_matrix_complex_memcpy(gsl_matrix_complex * dest, const gsl_matrix_complex * src);
0219 int gsl_matrix_complex_swap(gsl_matrix_complex * m1, gsl_matrix_complex * m2);
0220 int gsl_matrix_complex_tricpy(CBLAS_UPLO_t Uplo, CBLAS_DIAG_t Diag, gsl_matrix_complex * dest, const gsl_matrix_complex * src);
0221
0222 int gsl_matrix_complex_swap_rows(gsl_matrix_complex * m, const size_t i, const size_t j);
0223 int gsl_matrix_complex_swap_columns(gsl_matrix_complex * m, const size_t i, const size_t j);
0224 int gsl_matrix_complex_swap_rowcol(gsl_matrix_complex * m, const size_t i, const size_t j);
0225
0226 int gsl_matrix_complex_transpose (gsl_matrix_complex * m);
0227 int gsl_matrix_complex_transpose_memcpy (gsl_matrix_complex * dest, const gsl_matrix_complex * src);
0228 int gsl_matrix_complex_transpose_tricpy(CBLAS_UPLO_t Uplo_src, CBLAS_DIAG_t Diag, gsl_matrix_complex * dest, const gsl_matrix_complex * src);
0229
0230 int gsl_matrix_complex_conjtrans_memcpy (gsl_matrix_complex * dest, const gsl_matrix_complex * src);
0231
0232 int gsl_matrix_complex_equal (const gsl_matrix_complex * a, const gsl_matrix_complex * b);
0233
0234 int gsl_matrix_complex_isnull (const gsl_matrix_complex * m);
0235 int gsl_matrix_complex_ispos (const gsl_matrix_complex * m);
0236 int gsl_matrix_complex_isneg (const gsl_matrix_complex * m);
0237 int gsl_matrix_complex_isnonneg (const gsl_matrix_complex * m);
0238
0239 int gsl_matrix_complex_add (gsl_matrix_complex * a, const gsl_matrix_complex * b);
0240 int gsl_matrix_complex_sub (gsl_matrix_complex * a, const gsl_matrix_complex * b);
0241 int gsl_matrix_complex_mul_elements (gsl_matrix_complex * a, const gsl_matrix_complex * b);
0242 int gsl_matrix_complex_div_elements (gsl_matrix_complex * a, const gsl_matrix_complex * b);
0243 int gsl_matrix_complex_scale (gsl_matrix_complex * a, const gsl_complex x);
0244 int gsl_matrix_complex_scale_rows (gsl_matrix_complex * a, const gsl_vector_complex * x);
0245 int gsl_matrix_complex_scale_columns (gsl_matrix_complex * a, const gsl_vector_complex * x);
0246 int gsl_matrix_complex_add_constant (gsl_matrix_complex * a, const gsl_complex x);
0247 int gsl_matrix_complex_add_diagonal (gsl_matrix_complex * a, const gsl_complex x);
0248 int gsl_matrix_complex_conjugate (gsl_matrix_complex * a);
0249
0250
0251
0252
0253 int gsl_matrix_complex_get_row(gsl_vector_complex * v, const gsl_matrix_complex * m, const size_t i);
0254 int gsl_matrix_complex_get_col(gsl_vector_complex * v, const gsl_matrix_complex * m, const size_t j);
0255 int gsl_matrix_complex_set_row(gsl_matrix_complex * m, const size_t i, const gsl_vector_complex * v);
0256 int gsl_matrix_complex_set_col(gsl_matrix_complex * m, const size_t j, const gsl_vector_complex * v);
0257
0258
0259
0260
0261 INLINE_DECL gsl_complex gsl_matrix_complex_get(const gsl_matrix_complex * m, const size_t i, const size_t j);
0262 INLINE_DECL void gsl_matrix_complex_set(gsl_matrix_complex * m, const size_t i, const size_t j, const gsl_complex x);
0263
0264 INLINE_DECL gsl_complex * gsl_matrix_complex_ptr(gsl_matrix_complex * m, const size_t i, const size_t j);
0265 INLINE_DECL const gsl_complex * gsl_matrix_complex_const_ptr(const gsl_matrix_complex * m, const size_t i, const size_t j);
0266
0267 #ifdef HAVE_INLINE
0268
0269 INLINE_FUN
0270 gsl_complex
0271 gsl_matrix_complex_get(const gsl_matrix_complex * m,
0272 const size_t i, const size_t j)
0273 {
0274 #if GSL_RANGE_CHECK
0275 if (GSL_RANGE_COND(1))
0276 {
0277 gsl_complex zero = {{0,0}};
0278
0279 if (i >= m->size1)
0280 {
0281 GSL_ERROR_VAL("first index out of range", GSL_EINVAL, zero) ;
0282 }
0283 else if (j >= m->size2)
0284 {
0285 GSL_ERROR_VAL("second index out of range", GSL_EINVAL, zero) ;
0286 }
0287 }
0288 #endif
0289 return *(gsl_complex *)(m->data + 2*(i * m->tda + j)) ;
0290 }
0291
0292 INLINE_FUN
0293 void
0294 gsl_matrix_complex_set(gsl_matrix_complex * m,
0295 const size_t i, const size_t j, const gsl_complex x)
0296 {
0297 #if GSL_RANGE_CHECK
0298 if (GSL_RANGE_COND(1))
0299 {
0300 if (i >= m->size1)
0301 {
0302 GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
0303 }
0304 else if (j >= m->size2)
0305 {
0306 GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
0307 }
0308 }
0309 #endif
0310 *(gsl_complex *)(m->data + 2*(i * m->tda + j)) = x ;
0311 }
0312
0313 INLINE_FUN
0314 gsl_complex *
0315 gsl_matrix_complex_ptr(gsl_matrix_complex * m,
0316 const size_t i, const size_t j)
0317 {
0318 #if GSL_RANGE_CHECK
0319 if (GSL_RANGE_COND(1))
0320 {
0321 if (i >= m->size1)
0322 {
0323 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
0324 }
0325 else if (j >= m->size2)
0326 {
0327 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
0328 }
0329 }
0330 #endif
0331 return (gsl_complex *)(m->data + 2*(i * m->tda + j)) ;
0332 }
0333
0334 INLINE_FUN
0335 const gsl_complex *
0336 gsl_matrix_complex_const_ptr(const gsl_matrix_complex * m,
0337 const size_t i, const size_t j)
0338 {
0339 #if GSL_RANGE_CHECK
0340 if (GSL_RANGE_COND(1))
0341 {
0342 if (i >= m->size1)
0343 {
0344 GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
0345 }
0346 else if (j >= m->size2)
0347 {
0348 GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
0349 }
0350 }
0351 #endif
0352 return (const gsl_complex *)(m->data + 2*(i * m->tda + j)) ;
0353 }
0354
0355 #endif
0356
0357 __END_DECLS
0358
0359 #endif