1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
|
/* matrix/gsl_matrix_long_double.h
*
* Copyright (C) 1996, 1997, 1998, 1999, 2000, 2007 Gerard Jungman, Brian Gough
*
* This program is free software; you can redistribute it and/or modify
* it under the terms of the GNU General Public License as published by
* the Free Software Foundation; either version 3 of the License, or (at
* your option) any later version.
*
* This program is distributed in the hope that it will be useful, but
* WITHOUT ANY WARRANTY; without even the implied warranty of
* MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
* General Public License for more details.
*
* You should have received a copy of the GNU General Public License
* along with this program; if not, write to the Free Software
* Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 02110-1301, USA.
*/
#ifndef __GSL_MATRIX_LONG_DOUBLE_H__
#define __GSL_MATRIX_LONG_DOUBLE_H__
#include <stdlib.h>
#include <gsl/gsl_types.h>
#include <gsl/gsl_errno.h>
#include <gsl/gsl_inline.h>
#include <gsl/gsl_check_range.h>
#include <gsl/gsl_vector_long_double.h>
#undef __BEGIN_DECLS
#undef __END_DECLS
#ifdef __cplusplus
# define __BEGIN_DECLS extern "C" {
# define __END_DECLS }
#else
# define __BEGIN_DECLS /* empty */
# define __END_DECLS /* empty */
#endif
__BEGIN_DECLS
typedef struct
{
size_t size1;
size_t size2;
size_t tda;
long double * data;
gsl_block_long_double * block;
int owner;
} gsl_matrix_long_double;
typedef struct
{
gsl_matrix_long_double matrix;
} _gsl_matrix_long_double_view;
typedef _gsl_matrix_long_double_view gsl_matrix_long_double_view;
typedef struct
{
gsl_matrix_long_double matrix;
} _gsl_matrix_long_double_const_view;
typedef const _gsl_matrix_long_double_const_view gsl_matrix_long_double_const_view;
/* Allocation */
gsl_matrix_long_double *
gsl_matrix_long_double_alloc (const size_t n1, const size_t n2);
gsl_matrix_long_double *
gsl_matrix_long_double_calloc (const size_t n1, const size_t n2);
gsl_matrix_long_double *
gsl_matrix_long_double_alloc_from_block (gsl_block_long_double * b,
const size_t offset,
const size_t n1,
const size_t n2,
const size_t d2);
gsl_matrix_long_double *
gsl_matrix_long_double_alloc_from_matrix (gsl_matrix_long_double * m,
const size_t k1,
const size_t k2,
const size_t n1,
const size_t n2);
gsl_vector_long_double *
gsl_vector_long_double_alloc_row_from_matrix (gsl_matrix_long_double * m,
const size_t i);
gsl_vector_long_double *
gsl_vector_long_double_alloc_col_from_matrix (gsl_matrix_long_double * m,
const size_t j);
void gsl_matrix_long_double_free (gsl_matrix_long_double * m);
/* Views */
_gsl_matrix_long_double_view
gsl_matrix_long_double_submatrix (gsl_matrix_long_double * m,
const size_t i, const size_t j,
const size_t n1, const size_t n2);
_gsl_vector_long_double_view
gsl_matrix_long_double_row (gsl_matrix_long_double * m, const size_t i);
_gsl_vector_long_double_view
gsl_matrix_long_double_column (gsl_matrix_long_double * m, const size_t j);
_gsl_vector_long_double_view
gsl_matrix_long_double_diagonal (gsl_matrix_long_double * m);
_gsl_vector_long_double_view
gsl_matrix_long_double_subdiagonal (gsl_matrix_long_double * m, const size_t k);
_gsl_vector_long_double_view
gsl_matrix_long_double_superdiagonal (gsl_matrix_long_double * m, const size_t k);
_gsl_vector_long_double_view
gsl_matrix_long_double_subrow (gsl_matrix_long_double * m, const size_t i,
const size_t offset, const size_t n);
_gsl_vector_long_double_view
gsl_matrix_long_double_subcolumn (gsl_matrix_long_double * m, const size_t j,
const size_t offset, const size_t n);
_gsl_matrix_long_double_view
gsl_matrix_long_double_view_array (long double * base,
const size_t n1,
const size_t n2);
_gsl_matrix_long_double_view
gsl_matrix_long_double_view_array_with_tda (long double * base,
const size_t n1,
const size_t n2,
const size_t tda);
_gsl_matrix_long_double_view
gsl_matrix_long_double_view_vector (gsl_vector_long_double * v,
const size_t n1,
const size_t n2);
_gsl_matrix_long_double_view
gsl_matrix_long_double_view_vector_with_tda (gsl_vector_long_double * v,
const size_t n1,
const size_t n2,
const size_t tda);
_gsl_matrix_long_double_const_view
gsl_matrix_long_double_const_submatrix (const gsl_matrix_long_double * m,
const size_t i, const size_t j,
const size_t n1, const size_t n2);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_row (const gsl_matrix_long_double * m,
const size_t i);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_column (const gsl_matrix_long_double * m,
const size_t j);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_diagonal (const gsl_matrix_long_double * m);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_subdiagonal (const gsl_matrix_long_double * m,
const size_t k);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_superdiagonal (const gsl_matrix_long_double * m,
const size_t k);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_subrow (const gsl_matrix_long_double * m, const size_t i,
const size_t offset, const size_t n);
_gsl_vector_long_double_const_view
gsl_matrix_long_double_const_subcolumn (const gsl_matrix_long_double * m, const size_t j,
const size_t offset, const size_t n);
_gsl_matrix_long_double_const_view
gsl_matrix_long_double_const_view_array (const long double * base,
const size_t n1,
const size_t n2);
_gsl_matrix_long_double_const_view
gsl_matrix_long_double_const_view_array_with_tda (const long double * base,
const size_t n1,
const size_t n2,
const size_t tda);
_gsl_matrix_long_double_const_view
gsl_matrix_long_double_const_view_vector (const gsl_vector_long_double * v,
const size_t n1,
const size_t n2);
_gsl_matrix_long_double_const_view
gsl_matrix_long_double_const_view_vector_with_tda (const gsl_vector_long_double * v,
const size_t n1,
const size_t n2,
const size_t tda);
/* Operations */
void gsl_matrix_long_double_set_zero (gsl_matrix_long_double * m);
void gsl_matrix_long_double_set_identity (gsl_matrix_long_double * m);
void gsl_matrix_long_double_set_all (gsl_matrix_long_double * m, long double x);
int gsl_matrix_long_double_fread (FILE * stream, gsl_matrix_long_double * m) ;
int gsl_matrix_long_double_fwrite (FILE * stream, const gsl_matrix_long_double * m) ;
int gsl_matrix_long_double_fscanf (FILE * stream, gsl_matrix_long_double * m);
int gsl_matrix_long_double_fprintf (FILE * stream, const gsl_matrix_long_double * m, const char * format);
int gsl_matrix_long_double_memcpy(gsl_matrix_long_double * dest, const gsl_matrix_long_double * src);
int gsl_matrix_long_double_swap(gsl_matrix_long_double * m1, gsl_matrix_long_double * m2);
int gsl_matrix_long_double_swap_rows(gsl_matrix_long_double * m, const size_t i, const size_t j);
int gsl_matrix_long_double_swap_columns(gsl_matrix_long_double * m, const size_t i, const size_t j);
int gsl_matrix_long_double_swap_rowcol(gsl_matrix_long_double * m, const size_t i, const size_t j);
int gsl_matrix_long_double_transpose (gsl_matrix_long_double * m);
int gsl_matrix_long_double_transpose_memcpy (gsl_matrix_long_double * dest, const gsl_matrix_long_double * src);
long double gsl_matrix_long_double_max (const gsl_matrix_long_double * m);
long double gsl_matrix_long_double_min (const gsl_matrix_long_double * m);
void gsl_matrix_long_double_minmax (const gsl_matrix_long_double * m, long double * min_out, long double * max_out);
void gsl_matrix_long_double_max_index (const gsl_matrix_long_double * m, size_t * imax, size_t *jmax);
void gsl_matrix_long_double_min_index (const gsl_matrix_long_double * m, size_t * imin, size_t *jmin);
void gsl_matrix_long_double_minmax_index (const gsl_matrix_long_double * m, size_t * imin, size_t * jmin, size_t * imax, size_t * jmax);
int gsl_matrix_long_double_equal (const gsl_matrix_long_double * a, const gsl_matrix_long_double * b);
int gsl_matrix_long_double_isnull (const gsl_matrix_long_double * m);
int gsl_matrix_long_double_ispos (const gsl_matrix_long_double * m);
int gsl_matrix_long_double_isneg (const gsl_matrix_long_double * m);
int gsl_matrix_long_double_isnonneg (const gsl_matrix_long_double * m);
int gsl_matrix_long_double_add (gsl_matrix_long_double * a, const gsl_matrix_long_double * b);
int gsl_matrix_long_double_sub (gsl_matrix_long_double * a, const gsl_matrix_long_double * b);
int gsl_matrix_long_double_mul_elements (gsl_matrix_long_double * a, const gsl_matrix_long_double * b);
int gsl_matrix_long_double_div_elements (gsl_matrix_long_double * a, const gsl_matrix_long_double * b);
int gsl_matrix_long_double_scale (gsl_matrix_long_double * a, const double x);
int gsl_matrix_long_double_add_constant (gsl_matrix_long_double * a, const double x);
int gsl_matrix_long_double_add_diagonal (gsl_matrix_long_double * a, const double x);
/***********************************************************************/
/* The functions below are obsolete */
/***********************************************************************/
int gsl_matrix_long_double_get_row(gsl_vector_long_double * v, const gsl_matrix_long_double * m, const size_t i);
int gsl_matrix_long_double_get_col(gsl_vector_long_double * v, const gsl_matrix_long_double * m, const size_t j);
int gsl_matrix_long_double_set_row(gsl_matrix_long_double * m, const size_t i, const gsl_vector_long_double * v);
int gsl_matrix_long_double_set_col(gsl_matrix_long_double * m, const size_t j, const gsl_vector_long_double * v);
/***********************************************************************/
/* inline functions if you are using GCC */
INLINE_DECL long double gsl_matrix_long_double_get(const gsl_matrix_long_double * m, const size_t i, const size_t j);
INLINE_DECL void gsl_matrix_long_double_set(gsl_matrix_long_double * m, const size_t i, const size_t j, const long double x);
INLINE_DECL long double * gsl_matrix_long_double_ptr(gsl_matrix_long_double * m, const size_t i, const size_t j);
INLINE_DECL const long double * gsl_matrix_long_double_const_ptr(const gsl_matrix_long_double * m, const size_t i, const size_t j);
#ifdef HAVE_INLINE
INLINE_FUN
long double
gsl_matrix_long_double_get(const gsl_matrix_long_double * m, const size_t i, const size_t j)
{
#if GSL_RANGE_CHECK
if (GSL_RANGE_COND(1))
{
if (i >= m->size1)
{
GSL_ERROR_VAL("first index out of range", GSL_EINVAL, 0) ;
}
else if (j >= m->size2)
{
GSL_ERROR_VAL("second index out of range", GSL_EINVAL, 0) ;
}
}
#endif
return m->data[i * m->tda + j] ;
}
INLINE_FUN
void
gsl_matrix_long_double_set(gsl_matrix_long_double * m, const size_t i, const size_t j, const long double x)
{
#if GSL_RANGE_CHECK
if (GSL_RANGE_COND(1))
{
if (i >= m->size1)
{
GSL_ERROR_VOID("first index out of range", GSL_EINVAL) ;
}
else if (j >= m->size2)
{
GSL_ERROR_VOID("second index out of range", GSL_EINVAL) ;
}
}
#endif
m->data[i * m->tda + j] = x ;
}
INLINE_FUN
long double *
gsl_matrix_long_double_ptr(gsl_matrix_long_double * m, const size_t i, const size_t j)
{
#if GSL_RANGE_CHECK
if (GSL_RANGE_COND(1))
{
if (i >= m->size1)
{
GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
}
else if (j >= m->size2)
{
GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
}
}
#endif
return (long double *) (m->data + (i * m->tda + j)) ;
}
INLINE_FUN
const long double *
gsl_matrix_long_double_const_ptr(const gsl_matrix_long_double * m, const size_t i, const size_t j)
{
#if GSL_RANGE_CHECK
if (GSL_RANGE_COND(1))
{
if (i >= m->size1)
{
GSL_ERROR_NULL("first index out of range", GSL_EINVAL) ;
}
else if (j >= m->size2)
{
GSL_ERROR_NULL("second index out of range", GSL_EINVAL) ;
}
}
#endif
return (const long double *) (m->data + (i * m->tda + j)) ;
}
#endif
__END_DECLS
#endif /* __GSL_MATRIX_LONG_DOUBLE_H__ */
|