matrix.h
Go to the documentation of this file.
1 /*
2  Harmonic sinusoidal modelling and tools
3 
4  C++ code package for harmonic sinusoidal modelling and relevant signal processing.
5  Centre for Digital Music, Queen Mary, University of London.
6  This file copyright 2011 Wen Xue.
7 
8  This program is free software; you can redistribute it and/or
9  modify it under the terms of the GNU General Public License as
10  published by the Free Software Foundation; either version 2 of the
11  License, or (at your option) any later version.
12 */
13 #ifndef MatrixH
14 #define MatrixH
15 
22 #include "xcomplex.h"
23 #include "arrayalloc.h"
24 
25 //--Similar balance----------------------------------------------------------
26 void BalanceSim(int n, double** A);
27 
28 //--Choleski factorization---------------------------------------------------
29 int Choleski(int N, double** L, double** A);
30 
31 //--matrix copy--------------------------------------------------------------
32 double** Copy(int M, int N, double** Z, double** A, MList* List=0);
33  double** Copy(int M, int N, double** A, MList* List=0);
34  double** Copy(int N, double** Z, double ** A, MList* List=0);
35  double** Copy(int N, double ** A, MList* List=0);
36 cdouble** Copy(int M, int N, cdouble** Z, cdouble** A, MList* List=0);
37  cdouble** Copy(int M, int N, cdouble** A, MList* List=0);
38  cdouble** Copy(int N, cdouble** Z, cdouble** A, MList* List=0);
39  cdouble** Copy(int N, cdouble** A, MList* List=0);
40 double* Copy(int N, double* z, double* a, MList* List=0);
41  double* Copy(int N, double* a, MList* List=0);
42 cdouble* Copy(int N, cdouble* z, cdouble* a, MList* List=0);
43  cdouble* Copy(int N, cdouble* a, MList* List=0);
44 
45 //--matrix determinant calculation-------------------------------------------
46 double det(int N, double** A, int mode=0);
47 cdouble det(int N, cdouble** A, int mode=0);
48 
49 //--power methods for solving eigenproblems----------------------------------
50 int EigPower(int N, double& l, double* x, double** A, double ep=1e-6, int maxiter=50);
51 int EigPowerA(int N, double& l, double* x, double** A, double ep=1e-6, int maxiter=50);
52 int EigPowerI(int N, double& l, double* x, double** A, double ep=1e-6, int maxiter=50);
53 int EigPowerS(int N, double& l, double* x, double** A, double ep=1e-6, int maxiter=50);
54 int EigPowerWielandt(int N, double& m, double* u, double l, double* v, double** A, double ep=1e-06, int maxiter=50);
55 int EigenValues(int N, double** A, cdouble* ev);
56 int EigSym(int N, double** A, double* d, double** Q);
57 
58 //--Gaussian elimination for solving linear systems--------------------------
59 int GEB(int N, double* x, double** A, double* b);
60 int GESCP(int N, double* x, double** A, double*b);
61 void GExL(int N, double* x, double** L, double* a);
62 void GExLAdd(int N, double* x, double** L, double* a);
63 void GExL1(int N, double* x, double** L, double a);
64 void GExL1Add(int N, double* x, double** L, double a);
65 
66 /*
67  template GECP: Gaussian elimination with maximal column pivoting for solving linear system Ax=b
68 
69  In: matrix A[N][N], vector b[N]
70  Out: vector x[N]
71 
72  Returns 0 if successful. Contents of matrix A and vector b are destroyed on return.
73 */
74 template<class T, class Ta>int GECP(int N, T* x, Ta** A, T *b, Ta* logdet=0)
75 {
76  if (logdet) *logdet=1E-302; int c, p, ip, *rp=new int[N]; for (int i=0; i<N; i++) rp[i]=i;
77  Ta m;
78  //Gaussian eliminating
79  for (int i=0; i<N-1; i++)
80  {
81  p=i, ip=i+1;
82  while (ip<N){if (fabs(A[rp[ip]][i])>fabs(A[rp[p]][i])) p=ip; ip++;}
83  if (A[rp[p]][i]==0) {delete[] rp; return 1;}
84  if (p!=i) {c=rp[i]; rp[i]=rp[p]; rp[p]=c;}
85  for (int j=i+1; j<N; j++)
86  {
87  m=A[rp[j]][i]/A[rp[i]][i];
88  A[rp[j]][i]=0;
89  for (int k=i+1; k<N; k++) A[rp[j]][k]-=m*A[rp[i]][k];
90  b[rp[j]]-=m*b[rp[i]];
91  }
92  }
93  if (A[rp[N-1]][N-1]==0) {delete[] rp; return 1;}
94  //backward substitution
95  x[N-1]=b[rp[N-1]]/A[rp[N-1]][N-1];
96  for (int i=N-2; i>=0; i--)
97  {
98  x[i]=b[rp[i]]; for (int j=i+1; j<N; j++) x[i]-=A[rp[i]][j]*x[j]; x[i]/=A[rp[i]][i];
99  }
100  if (logdet){*logdet=log(fabs(A[rp[0]][0])); for (int n=1; n<N; n++) *logdet+=log(fabs(A[rp[n]][n]));}
101  delete[] rp;
102  return 0;
103 }//GECP
104 
105 //--inverse lower and upper triangular matrices------------------------------
106 double GILT(int N, double** A);
107 double GIUT(int N, double** A);
108 
109 //--inverse matrix calculation with gaussian elimination---------------------
110 double GICP(int N, double** X, double** A);
111 double GICP(int N, double** A);
112 cdouble GICP(int N, cdouble** A);
113 double GISCP(int N, double** X, double** A);
114 double GISCP(int N, double** A);
115 
116 //--Gaussian-Seidel method for solving linear systems------------------------
117 int GSI(int N, double* x0, double** A, double* b, double ep=1e-4, int maxiter=50);
118 
119 //Reduction to upper Heissenberg matrix by elimination with pivoting
120 void Hessenb(int n, double** A);
121 
122 //--Householder algorithm converting a matrix tridiagonal--------------------
123 void HouseHolder(int N, double** A);
124 void HouseHolder(int N, double** T, double** Q, double** A);
125 void HouseHolder(int n, double **A, double* d, double* sd);
126 
127 //--inner product------------------------------------------------------------
128 double Inner(int N, double* x, double* y);
129 cdouble Inner(int N, double* x, cdouble* y);
130 cdouble Inner(int N, cdouble* x, cdouble* y);
131 cdouble Inner(int N, cfloat* x, cdouble* y);
132 cfloat Inner(int N, cfloat* x, cfloat* y);
133 double Inner(int M, int N, double** X, double** Y);
134 
135 /*
136  template Inner: Inner product <xw, y>
137 
138  In: vectors x[N], w[N] and y[N]
139 
140  Returns inner product of xw and y.
141 */
142 template<class Tx, class Tw>cdouble Inner(int N, Tx* x, Tw* w, cdouble* y)
143 {
144  cdouble result=0;
145  for (int i=0; i<N; i++) result+=(x[i]*w[i])**y[i];
146  return result;
147 }//Inner
148 template<class Tx, class Tw>cdouble Inner(int N, Tx* x, Tw* w, double* y)
149 {
150  cdouble result=0;
151  for (int i=0; i<N; i++) result+=x[i]*w[i]*y[i];
152  return result;
153 }//Inner
154 
155 //--Jacobi iterative method for solving linear systems-----------------------
156 int JI(int N, double* x, double** A, double* b, double ep=1e-4, int maxiter=50);
157 
158 //--LDL factorization of a symmetric matrix----------------------------------
159 int LDL(int N, double** L, double* d, double** A);
160 
161 //--LQ factorization by Gram-Schmidt-----------------------------------------
162 void LQ_GS(int M, int N, double** A, double** L, double** Q);
163 
164 //--1-stage Least-square solution of overdetermined linear system------------
165 /*
166  template LSLinera: direct LS solution of A[M][N]x[N]=y[M], M>=N.
167 
168  In: matrix A[M][N], vector y[M]
169  Out: vector x[N]
170 
171  Returns the log determinant of AtA. Contents of matrix A and vector y are unchanged on return.
172 */
173 template <class T> T LSLinear(int M, int N, T* x, T** A, T* y)
174 {
175  MList* mlist=new MList;
176  T** AtA=MultiplyXtX(N, M, A, mlist);
177  T* Aty=MultiplyXty(N, M, A, y, mlist);
178  T result; GECP(N, x, AtA, Aty, &result);
179  delete mlist;
180  return result;
181 }//LSLinear
182 
183 //--2-stage Least-square solution of overdetermined linear system------------
184 void LSLinear2(int M, int N, double* x, double** A, double* y);
185 
186 //--LU factorization of a non-singular matrix--------------------------------
187 int LU_Direct(int mode, int N, double* diag, double** a);
188 int LU(int N, double** l, double** u, double** a); void LU_DiagL(int N, double** l, double* diagl, double** u, double** a);
189 int LU_PD(int N, double** b);
190 int LU_PD(int N, double** b, double* c);
191 double LUCP(double **a, int N, int *ind, int& parity, int mode=1);
192 
193 //--LU factorization method for solving a linear system----------------------
194 void LU(int N, double* x, double** A, double* y, int* ind=0);
195 
196 //--find maximum of vector---------------------------------------------------
197 int maxind(double* data, int from, int to);
198 
199 //--matrix linear combination------------------------------------------------
200 /*
201  template MultiAdd: matrix linear combination Z=X+aY
202 
203  In: matrices X[M][N], Y[M][N], scalar a
204  Out: matrix Z[M][N]
205 
206  Returns pointer to Z. Z is created anew if Z=0 is specified on start.
207 */
208 template<class Tz, class Tx, class Ty, class Ta> Tz** MultiAdd(int M, int N, Tz** Z, Tx** X, Ty** Y, Ta a, MList* List=0)
209 {
210  if (!Z){Allocate2(Tz, M, N, Z); if (List) List->Add(Z, 2);}
211  for (int i=0; i<M; i++) for (int j=0; j<N; j++) Z[i][j]=X[i][j]+a*Y[i][j];
212  return Z;
213 }//MultiAdd
214 
215 /*
216  template MultiAdd: vector linear combination z=x+ay
217 
218  In: vectors x[N], y[N], scalar a
219  Out: vector z[N]
220 
221  Returns pointer to z. z is created anew if z=0 is specified on start.
222 */
223 template<class Tz, class Tx, class Ty, class Ta> Tz* MultiAdd(int N, Tz* z, Tx* x, Ty* y, Ta a, MList* List=0)
224 {
225  if (!z){z=new Tz[N]; if (List) List->Add(z, 1);}
226  for (int i=0; i<N; i++) z[i]=x[i]+Ty(a)*y[i];
227  return z;
228 }//MultiAdd
229 template<class Tx, class Ty, class Ta> Tx* MultiAdd(int N, Tx* x, Ty* y, Ta a, MList* List=0)
230 {
231  return MultiAdd(N, (Tx*)0, x, y, a, List);
232 }//MultiAdd
233 
234 //--matrix multiplication by constant----------------------------------------
235 /*
236  template Multiply: matrix-constant multiplication Z=aX
237 
238  In: matrix X[M][N], scalar a
239  Out: matrix Z[M][N]
240 
241  Returns pointer to Z. Z is created anew if Z=0 is specified on start.
242 */
243 template<class Tz, class Tx, class Ta> Tz** Multiply(int M, int N, Tz** Z, Tx** X, Ta a, MList* List=0)
244 {
245  if (!Z){Allocate2(Tz, M, N, Z); if (List) List->Add(Z, 2);}
246  for (int i=0; i<M; i++) for (int j=0; j<N; j++) Z[i][j]=a*X[i][j];
247  return Z;
248 }//Multiply
249 
250 /*
251  template Multiply: vector-constant multiplication z=ax
252 
253  In: matrix x[N], scalar a
254  Out: matrix z[N]
255 
256  Returns pointer to z. z is created anew if z=0 is specified on start.
257 */
258 template<class Tz, class Tx, class Ta> Tz* Multiply(int N, Tz* z, Tx* x, Ta a, MList* List=0)
259 {
260  if (!z){z=new Tz[N]; if (List) List->Add(z, 1);}
261  for (int i=0; i<N; i++) z[i]=x[i]*a;
262  return z;
263 }//Multiply
264 template<class Tx, class Ta>Tx* Multiply(int N, Tx* x, Ta a, MList* List=0)
265 {
266  return Multiply(N, (Tx*)0, x, a, List);
267 }//Multiply
268 
269 //--matrix multiplication operations-----------------------------------------
270 /*
271  macro Multiply_full: matrix-matrix multiplication z=xy and multiplication-accumulation z=z+xy.
272 
273  Each expansion of the macro defines three function templates; two are named $MULTIPLY and do matrix
274  multiplication only; one is named $MULTIADD and accumulates the multiplicated result on top of a
275  specified matrix. One of the two $MULTIPLY functions allows using a pre-allocated matrix to accept
276  the matrix product, while the other directly allocates a new matrix and returns the pointer.
277 
278  Functions are named after their exact functions. For example, MultiplyXYc multiplies matrix X with
279  the complex conjugate of matrix Y, where postfix "c" attched to Y stands for conjugate. Likewise,
280  the postfix "t" stands for transpose, and "h" stnads for Hermitian (conjugate transpose).
281 
282  Three dimension arguments are needed by each function template. The first and last of the three are
283  the number of rows and columns, respectively, of the product (output) matrix. The middle one is the
284  "other", common dimension of both multiplier matrices.
285 */
286 #define Multiply_full(MULTIPLY, MULTIADD, xx, yy) \
287  template<class Tx, class Ty>Tx** MULTIPLY(int N, int M, int K, Tx** z, Tx** x, Ty** y, MList* List=0){ \
288  if (!z) {Allocate2(Tx, N, K, z); if (List) List->Add(z, 2);} \
289  for (int n=0; n<N; n++) for (int k=0; k<K; k++){ \
290  Tx zz=0; for (int m=0; m<M; m++) zz+=xx*yy; z[n][k]=zz;} return z;} \
291  template<class Tx, class Ty>Tx** MULTIPLY(int N, int M, int K, Tx** x, Ty** y, MList* List=0){ \
292  Tx** Allocate2(Tx, N, K, z); if (List) List->Add(z, 2); \
293  for (int n=0; n<N; n++) for (int k=0; k<K; k++){ \
294  Tx zz=0; for (int m=0; m<M; m++) zz+=xx*yy; z[n][k]=zz;} return z;} \
295  template<class Tx, class Ty>void MULTIADD(int N, int M, int K, Tx** z, Tx** x, Ty** y){ \
296  for (int n=0; n<N; n++) for (int k=0; k<K; k++){ \
297  Tx zz=0; for (int m=0; m<M; m++) zz+=xx*yy; z[n][k]+=zz;}}
298 Multiply_full(MultiplyXY, MultiAddXY, x[n][m], y[m][k]) //z[N][K]=x[N][M]y[M][K], identical z and x or y not allowed
299 Multiply_full(MultiplyXYc, MultiAddXYc, x[n][m], *y[m][k]) //z[N][K]=x[N][M](y*)[M][K], identical z and x or y not allowed
300 Multiply_full(MultiplyXYt, MultiAddXYt, x[n][m], y[k][m]) //z[N][K]=x[N][M]Yt[M][K], identical z and x or y not allowed
301 Multiply_full(MultiplyXYh, MultiAddXYh, x[n][m], *y[k][m]) //z[N][K]=x[N][M](Yt*)[M][K], identical z and x or y not allowed
302 
303 /*
304  macro Multiply_square: square matrix multiplication z=xy and multiplication-accumulation z=z+xy.
305  Identical z and x allowed for multiplication but not for multiplication-accumulation.
306 
307  Each expansion of the macro defines three function templates; two are named $MULTIPLY and do matrix
308  multiplication only; one is named $MULTIADD and accumulates the multiplicated result on top of a
309  specified matrix. One of the two $MULTIPLY functions allows using a pre-allocated matrix to accept
310  the matrix product, while the other directly allocates a new matrix and returns the pointer.
311 */
312 #define Multiply_square(MULTIPLY, MULTIADD) \
313  template<class T>T** MULTIPLY(int N, T** z, T** x, T** y, MList* List=0){ \
314  if (!z){Allocate2(T, N, N, z); if (List) List->Add(z, 2);} \
315  if (z!=x) MULTIPLY(N, N, N, z, x, y); else{ \
316  T* tmp=new T[N]; int sizeN=sizeof(T)*N; for (int i=0; i<N; i++){ \
317  MULTIPLY(1, N, N, &tmp, &x[i], y); memcpy(z[i], tmp, sizeN);} delete[] tmp;} \
318  return z;} \
319  template<class T>T** MULTIPLY(int N, T** x, T** y, MList* List=0){return MULTIPLY(N, N, N, x, y, List);}\
320  template<class T>void MULTIADD(int N, T** z, T** x, T** y){MULTIADD(N, N, N, z, x, y);}
321 Multiply_square(MultiplyXY, MultiAddXY)
322 
323 /*
324  macro Multiply_xx: matrix self multiplication z=xx and self-multiplication-accumulation z=z+xx.
325 
326  Each expansion of the macro defines three function templates; two are named $MULTIPLY and do matrix
327  multiplication only; one is named $MULTIADD and accumulates the multiplicated result on top of a
328  specified matrix. One of the two $MULTIPLY functions allows using a pre-allocated matrix to accept
329  the matrix product, while the other directly allocates a new matrix and returns the pointer.
330 */
331 #define Multiply_xx(MULTIPLY, MULTIADD, xx1, xx2, zzt) \
332  template<class T>T** MULTIPLY(int M, int N, T** z, T** x, MList* List=0){ \
333  if (!z){Allocate2(T, M, M, z); if (List) List->Add(z, 2);} \
334  for (int m=0; m<M; m++) for (int k=0; k<=m; k++){ \
335  T zz=0; for (int n=0; n<N; n++) zz+=xx1*xx2; z[m][k]=zz; if (m!=k) z[k][m]=zzt;} \
336  return z;} \
337  template<class T>T** MULTIPLY(int M, int N, T ** x, MList* List=0){ \
338  T** Allocate2(T, M, M, z); if (List) List->Add(z, 2); \
339  for (int m=0; m<M; m++) for (int k=0; k<=m; k++){ \
340  T zz=0; for (int n=0; n<N; n++) zz+=xx1*xx2; z[m][k]=zz; if (m!=k) z[k][m]=zzt;} \
341  return z;} \
342  template<class T>void MULTIADD(int M, int N, T** z, T** x){ \
343  for (int m=0; m<M; m++) for (int k=0; k<=m; k++){ \
344  T zz=0; for (int n=0; n<N; n++) zz+=xx1*xx2; z[m][k]+=zz; if (m!=k) z[k][m]+=zzt;}}
345 Multiply_xx(MultiplyXtX, MultiAddXtX, x[n][m], x[n][k], zz) //z[M][M]=xt[M][N]x[N][M], identical z and x NOT ALLOWED.
346 Multiply_xx(MultiplyXXt, MultiAddXXt, x[m][n], x[k][n], zz) //z[M][M]=X[M][N]xt[N][M], identical z and x NOT ALLOWED.
347 Multiply_xx(MultiplyXhX, MultiAddXhX, *x[n][m], x[n][k], *zz) //z[M][M]=(xt*)[M][N]x[N][M], identical z and x NOT ALLOWED.
348 Multiply_xx(MultiplyXXh, MultiAddXXh, x[m][n], *x[k][n], *zz) //z[M][M]=x[M][N](xt*)[N][M], identical z and x NOT ALLOWED.
349 Multiply_xx(MultiplyXcXt, MultiAddXcXt, *x[m][n], x[k][n], *zz) //z[M][M]=(x*)[M][N]xt[N][M], identical z and x NOT ALLOWED.
350 
351 //matrix-vector multiplication routines
352 double* MultiplyXy(int M, int N, double* z, double** x, double* y, MList* List=0);
353 double* MultiplyXy(int M, int N, double** x, double* y, MList* List=0);
354 cdouble* MultiplyXy(int M, int N, cdouble* z, cdouble** x, cdouble* y, MList* List=0);
355 cdouble* MultiplyXy(int M, int N, cdouble** x, cdouble* y, MList* List=0);
356 cdouble* MultiplyXy(int M, int N, cdouble* z, double** x, cdouble* y, MList* List=0);
357 cdouble* MultiplyXy(int M, int N, double** x, cdouble* y, MList* List=0);
358 double* MultiplyxY(int M, int N, double* zt, double* xt, double** y, MList* List=0);
359 double* MultiplyxY(int M, int N, double* xt, double** y, MList* List=0);
360 cdouble* MultiplyxY(int M, int N, cdouble* zt, cdouble* xt, cdouble** y, MList* List=0);
361 cdouble* MultiplyxY(int M, int N, cdouble* xt, cdouble** y, MList* List=0);
362 double* MultiplyXty(int M, int N, double* z, double** x, double* y, MList* List=0);
363 double* MultiplyXty(int M, int N, double** x, double* y, MList* List=0);
364 cdouble* MultiplyXty(int M, int N, cdouble* z, cdouble** x, cdouble* y, MList* List=0);
365 cdouble* MultiplyXty(int M, int N, cdouble** x, cdouble* y, MList* List=0);
366 cdouble* MultiplyXhy(int M, int N, cdouble* z, cdouble** x, cdouble* y, MList* List=0);
367 cdouble* MultiplyXhy(int M, int N, cdouble** x, cdouble* y, MList* List=0);
368 cdouble* MultiplyXcy(int M, int N, cdouble* z, cdouble** x, cfloat* y, MList* List=0);
369 cdouble* MultiplyXcy(int M, int N, cdouble** x, cfloat* y, MList* List=0);
370 cdouble* MultiplyXcy(int M, int N, cdouble* z, cdouble** x, cdouble* y, MList* List=0);
371 cdouble* MultiplyXcy(int M, int N, cdouble** x, cdouble* y, MList* List=0);
372 double* MultiplyxYt(int M, int N, double* zt, double* xt, double** y, MList* MList=0);
373 double* MultiplyxYt(int M, int N, double* xt, double** y, MList* MList=0);
374 
375 //--matrix norms-------------------------------------------------------------
376 double Norm1(int N, double** A);
377 
378 //--QR factorization---------------------------------------------------------
379 void QR_GS(int M, int N, double** A, double** Q, double** R);
380 void QR_householder(int M, int N, double** A, double** Q, double** R);
381 
382 //--QR factorization of a tridiagonal matrix---------------------------------
383 int QL(int N, double* d, double* sd, double** z);
384 int QR(int N, double **A, cdouble* ev);
385 
386 //--QU factorization of a matrix---------------------------------------------
387 void QU(int N, double** Q, double** A);
388 
389 //--Extract the real part----------------------------------------------------
390 double** Real(int M, int N, double** z, cdouble** x, MList* List);
391 double** Real(int M, int N, cdouble** x, MList* List);
392 
393 //--Finding roots of real polynomials----------------------------------------
394 int Roots(int N, double* p, double* rr, double* ri);
395 
396 //--Sor iterative method for solving linear systems--------------------------
397 int SorI(int N, double* x0, double** a, double* b, double w=1, double ep=1e-4, int maxiter=50);
398 
399 //--Submatrix----------------------------------------------------------------
400 void SetSubMatrix(double** z, double** x, int Y1, int Y, int X1, int X);
401 void SetSubMatrix(cdouble** z, cdouble** x, int Y1, int Y, int X1, int X);
402 cdouble** SubMatrix(cdouble** z, cdouble** x, int Y1, int Y, int X1, int X, MList* List=0);
403 cdouble** SubMatrix(cdouble** x, int Y1, int Y, int X1, int X, MList* List=0);
404 cdouble* SubVector(cdouble* z, cdouble* x, int X1, int X, MList* List=0);
405 cdouble* SubVector(cdouble* x, int X1, int X, MList* List=0);
406 
407 //--matrix transpose operation-----------------------------------------------
408 void transpose(int N, double** a);
409 void transpose(int N, cdouble** a);
410 double** transpose(int N, int M, double** ta, double** a, MList* List=0); //z[NM]=a[MN]'
411 double** transpose(int N, int M, double** a, MList* List=0);
412 
413 //--rotation matrix converting between vectors-------------------------------
414 double** Unitary(int N, double** P, double* x, double* y, MList* List=0);
415 double** Unitary(int N, double* x, double* y, MList* List=0);
416 cdouble** Unitary(int N, cdouble** P, cdouble* x, cdouble* y, MList* List=0);
417 cdouble** Unitary(int N, cdouble* x, cdouble* y, MList* List=0);
418 
419 #endif
void QR_householder(int M, int N, double **A, double **Q, double **R)
Definition: matrix.cpp:1950
void SetSubMatrix(double **z, double **x, int Y1, int Y, int X1, int X)
Definition: matrix.cpp:2104
cdouble * SubVector(cdouble *z, cdouble *x, int X1, int X, MList *List=0)
Definition: matrix.cpp:2142
void GExL1Add(int N, double *x, double **L, double a)
Definition: matrix.cpp:619
void GExLAdd(int N, double *x, double **L, double *a)
Definition: matrix.cpp:583
int EigPowerI(int N, double &l, double *x, double **A, double ep=1e-6, int maxiter=50)
Definition: matrix.cpp:329
int LU(int N, double **l, double **u, double **a)
Definition: matrix.cpp:1318
int EigPowerA(int N, double &l, double *x, double **A, double ep=1e-6, int maxiter=50)
Definition: matrix.cpp:295
void LSLinear2(int M, int N, double *x, double **A, double *y)
Definition: matrix.cpp:1277
double Inner(int N, double *x, double *y)
Definition: matrix.cpp:1120
Definition: arrayalloc.h:83
void QU(int N, double **Q, double **A)
Definition: matrix.cpp:1988
double GILT(int N, double **A)
Definition: matrix.cpp:730
int Choleski(int N, double **L, double **A)
Definition: matrix.cpp:89
int GSI(int N, double *x0, double **A, double *b, double ep=1e-4, int maxiter=50)
Definition: matrix.cpp:860
double det(int N, double **A, int mode=0)
Definition: matrix.cpp:179
double ** Copy(int M, int N, double **Z, double **A, MList *List=0)
Definition: matrix.cpp:120
double GIUT(int N, double **A)
Definition: matrix.cpp:756
int EigPowerWielandt(int N, double &m, double *u, double l, double *v, double **A, double ep=1e-06, int maxiter=50)
Definition: matrix.cpp:403
double ** Real(int M, int N, double **z, cdouble **x, MList *List)
Definition: matrix.cpp:2018
double LUCP(double **a, int N, int *ind, int &parity, int mode=1)
Definition: matrix.cpp:1555
void GExL1(int N, double *x, double **L, double a)
Definition: matrix.cpp:599
int SorI(int N, double *x0, double **a, double *b, double w=1, double ep=1e-4, int maxiter=50)
Definition: matrix.cpp:2070
void BalanceSim(int n, double **A)
Definition: matrix.cpp:31
int EigPowerS(int N, double &l, double *x, double **A, double ep=1e-6, int maxiter=50)
Definition: matrix.cpp:369
int JI(int N, double *x, double **A, double *b, double ep=1e-4, int maxiter=50)
Definition: matrix.cpp:1179
int EigenValues(int N, double **A, cdouble *ev)
Definition: matrix.cpp:442
int QL(int N, double *d, double *sd, double **z)
Definition: matrix.cpp:1723
void HouseHolder(int N, double **A)
Definition: matrix.cpp:948
int GESCP(int N, double *x, double **A, double *b)
Definition: matrix.cpp:518
double ** Unitary(int N, double **P, double *x, double *y, MList *List=0)
Definition: matrix.cpp:2205
int LU_PD(int N, double **b)
Definition: matrix.cpp:1458
void MultiAdd(double *out, double *in, double *adder, int count, double mul)
Definition: procedures.cpp:2148
int maxind(double *data, int from, int to)
Definition: matrix.cpp:1640
cdouble ** SubMatrix(cdouble **z, cdouble **x, int Y1, int Y, int X1, int X, MList *List=0)
Definition: matrix.cpp:2122
void Hessenb(int n, double **A)
Definition: matrix.cpp:892
Definition: xcomplex.h:26
int EigSym(int N, double **A, double *d, double **Q)
Definition: matrix.cpp:457
int EigPower(int N, double &l, double *x, double **A, double ep=1e-6, int maxiter=50)
Definition: matrix.cpp:262
int QR(int N, double **A, cdouble *ev)
Definition: matrix.cpp:1791
double GISCP(int N, double **X, double **A)
Definition: matrix.cpp:841
int GEB(int N, double *x, double **A, double *b)
Definition: matrix.cpp:476
int LDL(int N, double **L, double *d, double **A)
Definition: matrix.cpp:1216
void transpose(int N, double **a)
Definition: matrix.cpp:2163
void LQ_GS(int M, int N, double **A, double **L, double **Q)
Definition: matrix.cpp:1247
double GICP(int N, double **A)
Definition: matrix.cpp:636
int LU_Direct(int mode, int N, double *diag, double **a)
Definition: matrix.cpp:1403
void QR_GS(int M, int N, double **A, double **Q, double **R)
Definition: matrix.cpp:1922
void GExL(int N, double *x, double **L, double *a)
Definition: matrix.cpp:565