-
Notifications
You must be signed in to change notification settings - Fork 142
/
acb_mat.h
479 lines (348 loc) · 14.3 KB
/
acb_mat.h
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
351
352
353
354
355
356
357
358
359
360
361
362
363
364
365
366
367
368
369
370
371
372
373
374
375
376
377
378
379
380
381
382
383
384
385
386
387
388
389
390
391
392
393
394
395
396
397
398
399
400
401
402
403
404
405
406
407
408
409
410
411
412
413
414
415
416
417
418
419
420
421
422
423
424
425
426
427
428
429
430
431
432
433
434
435
436
437
438
439
440
441
442
443
444
445
446
447
448
449
450
451
452
453
454
455
456
457
458
459
460
461
462
463
464
465
466
467
468
469
470
471
472
473
474
475
476
477
478
/*
Copyright (C) 2012 Fredrik Johansson
This file is part of Arb.
Arb is free software: you can redistribute it and/or modify it under
the terms of the GNU Lesser General Public License (LGPL) as published
by the Free Software Foundation; either version 2.1 of the License, or
(at your option) any later version. See <http://www.gnu.org/licenses/>.
*/
#ifndef ACB_MAT_H
#define ACB_MAT_H
#ifdef ACB_MAT_INLINES_C
#define ACB_MAT_INLINE
#else
#define ACB_MAT_INLINE static __inline__
#endif
#include <stdio.h>
#include "flint/fmpz_mat.h"
#include "flint/fmpq_mat.h"
#include "arb.h"
#include "acb.h"
#include "arb_mat.h"
#include "acb_poly.h"
#ifdef __cplusplus
extern "C" {
#endif
typedef struct
{
acb_ptr entries;
slong r;
slong c;
acb_ptr * rows;
}
acb_mat_struct;
typedef acb_mat_struct acb_mat_t[1];
#define acb_mat_entry(mat,i,j) ((mat)->rows[i] + (j))
#define acb_mat_nrows(mat) ((mat)->r)
#define acb_mat_ncols(mat) ((mat)->c)
ACB_MAT_INLINE acb_ptr
acb_mat_entry_ptr(acb_mat_t mat, slong i, slong j)
{
return acb_mat_entry(mat, i, j);
}
/* Memory management */
void acb_mat_init(acb_mat_t mat, slong r, slong c);
void acb_mat_clear(acb_mat_t mat);
ACB_MAT_INLINE void
acb_mat_swap(acb_mat_t mat1, acb_mat_t mat2)
{
acb_mat_struct t = *mat1;
*mat1 = *mat2;
*mat2 = t;
}
ACB_MAT_INLINE void
acb_mat_swap_entrywise(acb_mat_t mat1, acb_mat_t mat2)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(mat1); i++)
for (j = 0; j < acb_mat_ncols(mat1); j++)
acb_swap(acb_mat_entry(mat2, i, j), acb_mat_entry(mat1, i, j));
}
/* Window matrices */
void acb_mat_window_init(acb_mat_t window, const acb_mat_t mat, slong r1, slong c1, slong r2, slong c2);
ARB_MAT_INLINE void
acb_mat_window_clear(acb_mat_t window)
{
flint_free(window->rows);
}
/* Conversions */
void acb_mat_set(acb_mat_t dest, const acb_mat_t src);
void acb_mat_set_fmpz_mat(acb_mat_t dest, const fmpz_mat_t src);
void acb_mat_set_round_fmpz_mat(acb_mat_t dest, const fmpz_mat_t src, slong prec);
void acb_mat_set_fmpq_mat(acb_mat_t dest, const fmpq_mat_t src, slong prec);
void acb_mat_set_arb_mat(acb_mat_t dest, const arb_mat_t src);
void acb_mat_set_round_arb_mat(acb_mat_t dest, const arb_mat_t src, slong prec);
/* Random generation */
void acb_mat_randtest(acb_mat_t mat, flint_rand_t state, slong prec, slong mag_bits);
void acb_mat_randtest_eig(acb_mat_t A, flint_rand_t state, acb_srcptr E, slong prec);
/* I/O */
void acb_mat_fprintd(FILE * file, const acb_mat_t mat, slong digits);
ACB_MAT_INLINE void
acb_mat_printd(const acb_mat_t mat, slong digits)
{
acb_mat_fprintd(stdout, mat, digits);
}
/* Comparisons */
int acb_mat_eq(const acb_mat_t mat1, const acb_mat_t mat2);
int acb_mat_ne(const acb_mat_t mat1, const acb_mat_t mat2);
int acb_mat_equal(const acb_mat_t mat1, const acb_mat_t mat2);
int acb_mat_overlaps(const acb_mat_t mat1, const acb_mat_t mat2);
int acb_mat_contains(const acb_mat_t mat1, const acb_mat_t mat2);
int acb_mat_contains_fmpq_mat(const acb_mat_t mat1, const fmpq_mat_t mat2);
int acb_mat_contains_fmpz_mat(const acb_mat_t mat1, const fmpz_mat_t mat2);
int acb_mat_is_real(const acb_mat_t mat);
ACB_MAT_INLINE int
acb_mat_is_empty(const acb_mat_t mat)
{
return (mat->r == 0) || (mat->c == 0);
}
ACB_MAT_INLINE int
acb_mat_is_square(const acb_mat_t mat)
{
return (mat->r == mat->c);
}
int acb_mat_is_exact(const acb_mat_t mat);
int acb_mat_is_zero(const acb_mat_t mat);
int acb_mat_is_finite(const acb_mat_t mat);
int acb_mat_is_triu(const acb_mat_t mat);
int acb_mat_is_tril(const acb_mat_t mat);
ACB_MAT_INLINE int
acb_mat_is_diag(const acb_mat_t mat)
{
return acb_mat_is_tril(mat) && acb_mat_is_triu(mat);
}
/* Radius and interval operations */
ACB_MAT_INLINE void
acb_mat_get_mid(acb_mat_t B, const acb_mat_t A)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
{
for (j = 0; j < acb_mat_ncols(A); j++)
{
arb_get_mid_arb(acb_realref(acb_mat_entry(B, i, j)), acb_realref(acb_mat_entry(A, i, j)));
arb_get_mid_arb(acb_imagref(acb_mat_entry(B, i, j)), acb_imagref(acb_mat_entry(A, i, j)));
}
}
}
ACB_MAT_INLINE void
acb_mat_add_error_mag(acb_mat_t mat, const mag_t err)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(mat); i++)
for (j = 0; j < acb_mat_ncols(mat); j++)
acb_add_error_mag(acb_mat_entry(mat, i, j), err);
}
/* Special matrices */
void acb_mat_zero(acb_mat_t mat);
void acb_mat_one(acb_mat_t mat);
void acb_mat_ones(acb_mat_t mat);
void acb_mat_indeterminate(acb_mat_t mat);
void acb_mat_dft(acb_mat_t res, int kind, slong prec);
void acb_mat_transpose(acb_mat_t mat1, const acb_mat_t mat2);
void acb_mat_conjugate(acb_mat_t mat1, const acb_mat_t mat2);
ACB_MAT_INLINE void
acb_mat_conjugate_transpose(acb_mat_t mat1, const acb_mat_t mat2)
{
acb_mat_transpose(mat1, mat2);
acb_mat_conjugate(mat1, mat1);
}
/* Norms */
void acb_mat_bound_inf_norm(mag_t b, const acb_mat_t A);
void acb_mat_frobenius_norm(arb_t res, const acb_mat_t A, slong prec);
void acb_mat_bound_frobenius_norm(mag_t b, const acb_mat_t A);
/* Arithmetic */
void acb_mat_neg(acb_mat_t dest, const acb_mat_t src);
void acb_mat_add(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_sub(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_mul_classical(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_mul_threaded(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_mul_reorder(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_mul(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_mul_entrywise(acb_mat_t res, const acb_mat_t mat1, const acb_mat_t mat2, slong prec);
void acb_mat_sqr_classical(acb_mat_t res, const acb_mat_t mat, slong prec);
void acb_mat_sqr(acb_mat_t res, const acb_mat_t mat, slong prec);
void acb_mat_pow_ui(acb_mat_t B, const acb_mat_t A, ulong exp, slong prec);
/* Scalar arithmetic */
ACB_MAT_INLINE void
acb_mat_scalar_mul_2exp_si(acb_mat_t B, const acb_mat_t A, slong c)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_mul_2exp_si(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c);
}
ACB_MAT_INLINE void
acb_mat_scalar_addmul_si(acb_mat_t B, const acb_mat_t A, slong c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_addmul_si(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_mul_si(acb_mat_t B, const acb_mat_t A, slong c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_mul_si(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_div_si(acb_mat_t B, const acb_mat_t A, slong c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_div_si(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_addmul_fmpz(acb_mat_t B, const acb_mat_t A, const fmpz_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_addmul_fmpz(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_mul_fmpz(acb_mat_t B, const acb_mat_t A, const fmpz_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_mul_fmpz(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_div_fmpz(acb_mat_t B, const acb_mat_t A, const fmpz_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_div_fmpz(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_addmul_acb(acb_mat_t B, const acb_mat_t A, const acb_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_addmul(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_mul_acb(acb_mat_t B, const acb_mat_t A, const acb_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_mul(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_div_acb(acb_mat_t B, const acb_mat_t A, const acb_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_div(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_addmul_arb(acb_mat_t B, const acb_mat_t A, const arb_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_addmul_arb(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_mul_arb(acb_mat_t B, const acb_mat_t A, const arb_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_mul_arb(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
ACB_MAT_INLINE void
acb_mat_scalar_div_arb(acb_mat_t B, const acb_mat_t A, const arb_t c, slong prec)
{
slong i, j;
for (i = 0; i < acb_mat_nrows(A); i++)
for (j = 0; j < acb_mat_ncols(A); j++)
acb_div_arb(acb_mat_entry(B, i, j), acb_mat_entry(A, i, j), c, prec);
}
/* Solving */
ACB_MAT_INLINE void
acb_mat_swap_rows(acb_mat_t mat, slong * perm, slong r, slong s)
{
if (r != s)
{
acb_ptr u;
slong t;
if (perm != NULL)
{
t = perm[s];
perm[s] = perm[r];
perm[r] = t;
}
u = mat->rows[s];
mat->rows[s] = mat->rows[r];
mat->rows[r] = u;
}
}
slong acb_mat_find_pivot_partial(const acb_mat_t mat,
slong start_row, slong end_row, slong c);
void acb_mat_solve_tril_classical(acb_mat_t X, const acb_mat_t L, const acb_mat_t B, int unit, slong prec);
void acb_mat_solve_tril_recursive(acb_mat_t X, const acb_mat_t L, const acb_mat_t B, int unit, slong prec);
void acb_mat_solve_tril(acb_mat_t X, const acb_mat_t L, const acb_mat_t B, int unit, slong prec);
void acb_mat_solve_triu_classical(acb_mat_t X, const acb_mat_t U, const acb_mat_t B, int unit, slong prec);
void acb_mat_solve_triu_recursive(acb_mat_t X, const acb_mat_t U, const acb_mat_t B, int unit, slong prec);
void acb_mat_solve_triu(acb_mat_t X, const acb_mat_t U, const acb_mat_t B, int unit, slong prec);
int acb_mat_lu_classical(slong * P, acb_mat_t LU, const acb_mat_t A, slong prec);
int acb_mat_lu_recursive(slong * P, acb_mat_t LU, const acb_mat_t A, slong prec);
int acb_mat_lu(slong * P, acb_mat_t LU, const acb_mat_t A, slong prec);
void acb_mat_solve_lu_precomp(acb_mat_t X, const slong * perm,
const acb_mat_t A, const acb_mat_t B, slong prec);
int acb_mat_solve_lu(acb_mat_t X, const acb_mat_t A, const acb_mat_t B, slong prec);
int acb_mat_solve(acb_mat_t X, const acb_mat_t A, const acb_mat_t B, slong prec);
int acb_mat_solve_precond(acb_mat_t X, const acb_mat_t A, const acb_mat_t B, slong prec);
void acb_mat_approx_mul(acb_mat_t C, const acb_mat_t A, const acb_mat_t B, slong prec);
void acb_mat_approx_solve_triu(acb_mat_t X, const acb_mat_t U, const acb_mat_t B, int unit, slong prec);
void acb_mat_approx_solve_tril(acb_mat_t X, const acb_mat_t L, const acb_mat_t B, int unit, slong prec);
int acb_mat_approx_lu(slong * P, acb_mat_t LU, const acb_mat_t A, slong prec);
void acb_mat_approx_solve_lu_precomp(acb_mat_t X, const slong * perm, const acb_mat_t A, const acb_mat_t B, slong prec);
int acb_mat_approx_solve(acb_mat_t X, const acb_mat_t A, const acb_mat_t B, slong prec);
int acb_mat_approx_inv(acb_mat_t X, const acb_mat_t A, slong prec);
int acb_mat_inv(acb_mat_t X, const acb_mat_t A, slong prec);
void acb_mat_det_lu(acb_t det, const acb_mat_t A, slong prec);
void acb_mat_det_precond(acb_t det, const acb_mat_t A, slong prec);
void acb_mat_det(acb_t det, const acb_mat_t A, slong prec);
/* Eigenvalues and eigenvectors */
int acb_mat_approx_eig_qr(acb_ptr E, acb_mat_t L, acb_mat_t R, const acb_mat_t A, const mag_t tol, slong maxiter, slong prec);
void acb_mat_eig_global_enclosure(mag_t eps, const acb_mat_t A, acb_srcptr E, const acb_mat_t R, slong prec);
void acb_mat_eig_enclosure_rump(acb_t lambda, acb_mat_t J, acb_mat_t X, const acb_mat_t A,
const acb_t lambda_approx, const acb_mat_t X_approx, slong prec);
int acb_mat_eig_simple_rump(acb_ptr E, acb_mat_t L, acb_mat_t R,
const acb_mat_t A, acb_srcptr E_approx, const acb_mat_t R_approx, slong prec);
int acb_mat_eig_simple_vdhoeven_mourrain(acb_ptr E, acb_mat_t L, acb_mat_t R,
const acb_mat_t A, acb_srcptr E_approx, const acb_mat_t R_approx, slong prec);
int acb_mat_eig_simple(acb_ptr E, acb_mat_t L, acb_mat_t R,
const acb_mat_t A, acb_srcptr E_approx, const acb_mat_t R_approx, slong prec);
int acb_mat_eig_multiple_rump(acb_ptr E, const acb_mat_t A, acb_srcptr E_approx, const acb_mat_t R_approx, slong prec);
int acb_mat_eig_multiple(acb_ptr E, const acb_mat_t A, acb_srcptr E_approx, const acb_mat_t R_approx, slong prec);
/* Special functions */
void acb_mat_exp_taylor_sum(acb_mat_t S, const acb_mat_t A, slong N, slong prec);
void acb_mat_exp(acb_mat_t B, const acb_mat_t A, slong prec);
void _acb_mat_charpoly(acb_ptr poly, const acb_mat_t mat, slong prec);
void acb_mat_charpoly(acb_poly_t poly, const acb_mat_t mat, slong prec);
void _acb_mat_companion(acb_mat_t mat, acb_srcptr poly, slong prec);
void acb_mat_companion(acb_mat_t mat, const acb_poly_t poly, slong prec);
void acb_mat_trace(acb_t trace, const acb_mat_t mat, slong prec);
void _acb_mat_diag_prod(acb_t res, const acb_mat_t A, slong a, slong b, slong prec);
void acb_mat_diag_prod(acb_t res, const acb_mat_t A, slong prec);
ACB_MAT_INLINE slong
acb_mat_allocated_bytes(const acb_mat_t x)
{
return _acb_vec_allocated_bytes(x->entries, x->r * x->c) + x->r * sizeof(acb_ptr);
}
#ifdef __cplusplus
}
#endif
#endif