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