1 /***************************************************************************
2  *   Copyright (C) 2009 by BUI Quang Minh   *
3  *   minh.bui@univie.ac.at   *
4  *                                                                         *
5  *   This program is free software; you can redistribute it and/or modify  *
6  *   it under the terms of the GNU General Public License as published by  *
7  *   the Free Software Foundation; either version 2 of the License, or     *
8  *   (at your option) any later version.                                   *
9  *                                                                         *
10  *   This program is distributed in the hope that it will be useful,       *
11  *   but WITHOUT ANY WARRANTY; without even the implied warranty of        *
12  *   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the         *
13  *   GNU General Public License for more details.                          *
14  *                                                                         *
15  *   You should have received a copy of the GNU General Public License     *
16  *   along with this program; if not, write to the                         *
17  *   Free Software Foundation, Inc.,                                       *
18  *   59 Temple Place - Suite 330, Boston, MA  02111-1307, USA.             *
19  ***************************************************************************/
20 #include "eigendecomposition.h"
21 #include "optimization.h"
22 #include <cmath>
23 #include <string.h>
24 #include <iostream>
25 #include <stdlib.h>
26 #include "tools.h"
27 
28 using namespace std;
29 
EigenDecomposition()30 EigenDecomposition::EigenDecomposition()
31 {
32 	total_num_subst = 1.0;
33 	normalize_matrix = true;
34     ignore_state_freq = false;
35 }
36 
eigensystem(double ** rate_params,double * state_freq,double * eval,double ** evec,double ** inv_evec,int num_state)37 void EigenDecomposition::eigensystem(
38 	double **rate_params, double *state_freq,
39 	double *eval, double **evec, double **inv_evec, int num_state)
40 {
41 	double *forg = new double[num_state];
42 	double *evali = new double[num_state];
43 	double *new_forg = new double[num_state];
44 	double *eval_new = new double[num_state];
45 	double **a = (double**)new double[num_state];
46 	double **b = (double**)new double[num_state];
47 	double **evec_new = (double**)new double[num_state];
48 	int *ordr = new int[num_state + 1];
49 	int i, j, k, error, new_num, inew, jnew;
50 	double zero;
51 
52 
53 	for (i=0; i < num_state; i++)
54 		a[i] = new double[num_state];
55 	for (i=0; i < num_state; i++)
56 		b[i] = new double[num_state];
57 	for (i=0; i < num_state; i++)
58 		evec_new[i] = new double[num_state];
59 
60 	/* get relative transition matrix and frequencies */
61 	memcpy(forg, state_freq, num_state * sizeof(double));
62     // BQM 2015-09-07: normalize state frequencies to 1
63     double sum = 0.0;
64     for (i = 0; i < num_state; i++)
65         sum += forg[i];
66     sum = 1.0/sum;
67     for (i = 0; i < num_state; i++)
68         forg[i] *= sum;
69 
70 	for (i = 0; i < num_state; i++)
71 		memcpy(a[i], rate_params[i], num_state * sizeof(double));
72 
73 	//rtfdata(a, forg, num_state);
74 	//    write (a, forg);
75 
76 	computeRateMatrix(a, forg, num_state); /* make 1 PAM rate matrix */
77 
78 	/* copy a to b */
79 	for (i = 0; i < num_state; i++)
80 		for (j = 0; j < num_state; j++)
81 			b[i][j] = a[i][j];
82 
83 	eliminateZero(b, forg, num_state, a, new_forg, new_num);
84 
85 	elmhes(a, ordr, new_num); /* compute eigenvalues and eigenvectors */
86 	//    writeInt (ordr);
87 
88 	eltran(a, evec_new, ordr, new_num);
89 
90 	//  writeMat (evec);
91 
92 	hqr2(new_num, 1, new_num, a, evec_new, eval_new, evali);
93 
94 
95 	// now get back eigen
96 	//for (i = 0,inew = 0; i < num_state; i++)
97 	for (i = num_state-1,inew = new_num-1; i >= 0; i--)
98 		eval[i] = (forg[i] > ZERO_FREQ) ? eval_new[inew--] : 0;
99 
100 	// calculate the actual eigenvectors of Q and its inverse matrix
101 	//for (i = 0, inew = 0; i < num_state; i++)
102 	for (i = num_state-1,inew = new_num-1; i >= 0; i--)
103 		if (forg[i] > ZERO_FREQ) {
104 			for (j = num_state-1, jnew = new_num-1; j >= 0; j--)
105 				if (forg[j] > ZERO_FREQ) {
106 					evec[i][j] = evec_new[inew][jnew];
107 					jnew--;
108 				} else {
109 					evec[i][j] = (i == j);
110 				}
111  			inew--;
112 		} else
113 		for (j=0; j < num_state; j++) {
114 			evec[i][j] = (i==j);
115 		}
116 
117 	/* check eigenvalue equation */
118 	error = 0;
119 	for (j = 0; j < num_state; j++) {
120 		for (i = 0, zero = 0.0; i < num_state; i++) {
121 			for (k = 0; k < num_state; k++) zero += b[i][k] * evec[k][j];
122 			zero -= eval[j] * evec[i][j];
123 			if (fabs(zero) > 1.0e-5) {
124 				error = 1;
125 				break;
126 			}
127 		}
128 	}
129 	if (error) {
130 		cout << "\nWARNING: Eigensystem doesn't satisfy eigenvalue equation!\n";
131 		cout << "Rate matrix R: " << endl;
132 		for (i = 0; i < num_state; i++) {
133 			for (j = 0; j < num_state; j++) cout << rate_params[i][j] << " ";
134 			cout << endl;
135 		}
136 		cout << "State frequencies: " << endl;
137 		for (i = 0; i < num_state; i++) cout << state_freq[i] << " ";
138 		cout << endl;
139 	}
140 
141 	for (i=num_state-1; i>= 0; i--)
142 		delete [] evec_new[i];
143 	for (i=num_state-1; i>= 0; i--)
144 		delete [] b[i];
145 	for (i=num_state-1; i>= 0; i--)
146 		delete [] a[i];
147 	delete [] ordr;
148 	delete [] evec_new;
149 	delete [] b;
150 	delete [] a;
151 	delete [] eval_new;
152 	delete [] new_forg;
153 	delete [] evali;
154 	delete [] forg;
155 
156 	luinverse(evec, inv_evec, num_state); /* inverse eigenvectors are in Ievc */
157 	//checkevector(evec, inv_evec, num_state); /* check whether inversion was OK */
158 
159 } /* eigensystem */
160 
161 
eigensystem_sym(double ** rate_params,double * state_freq,double * eval,double * evec,double * inv_evec,int num_state)162 void EigenDecomposition::eigensystem_sym(double **rate_params, double *state_freq,
163 	double *eval, double *evec, double *inv_evec, int num_state)
164 {
165 	double *forg = new double[num_state];
166 	double *new_forg = new double[num_state];
167 	double *forg_sqrt = new double[num_state];
168 	double *off_diag = new double[num_state];
169 	double *eval_new = new double[num_state];
170 	double **a = (double**)new double[num_state];
171 	double **b = (double**)new double[num_state];
172 	int i, j, k, new_num, inew, jnew;
173 	double error = 0.0;
174 	double zero;
175 
176 	for (i=0; i < num_state; i++)
177 		a[i] = new double[num_state];
178 	for (i=0; i < num_state; i++)
179 		b[i] = new double[num_state];
180 
181 	/* get relative transition matrix and frequencies */
182 	memcpy(forg, state_freq, num_state * sizeof(double));
183 
184     // BQM 2015-09-07: normalize state frequencies to 1
185     double sum = 0.0;
186     for (i = 0; i < num_state; i++)
187         sum += forg[i];
188     sum = 1.0/sum;
189     for (i = 0; i < num_state; i++)
190         forg[i] *= sum;
191 
192 
193 	for (i = 0; i < num_state; i++)
194 		memcpy(a[i], rate_params[i], num_state * sizeof(double));
195 
196 	//rtfdata(a, forg, num_state);
197 	//    write (a, forg);
198 
199 	computeRateMatrix(a, forg, num_state); /* make 1 PAM rate matrix */
200 
201 	/* copy a to b */
202 	for (i = 0; i < num_state; i++)
203 		for (j = 0; j < num_state; j++)
204 			b[i][j] = a[i][j];
205 
206 	eliminateZero(b, forg, num_state, a, new_forg, new_num);
207 
208 	symmetrizeRateMatrix(a, new_forg, forg_sqrt, new_num);
209 
210 	// make this matrix tridiagonal
211 	tred2(a, new_num, eval_new, off_diag);
212 	// compute eigenvalues and eigenvectors
213 	tqli(eval_new, off_diag, new_num, a);
214 
215     // make sure that all eval are non-positive
216 
217     for (i = 0; i < new_num; i++)
218         ASSERT(eval_new[i] <= 0.01);
219 
220 	// now get back eigen
221 	//for (i = 0,inew = 0; i < num_state; i++)
222 	for (i = num_state-1,inew = new_num-1; i >= 0; i--)
223 		eval[i] = (forg[i] > ZERO_FREQ) ? eval_new[inew--] : 0;
224 
225 	ASSERT(inew == -1);
226 	// calculate the actual eigenvectors of Q and its inverse matrix
227 	//for (i = 0, inew = 0; i < num_state; i++)
228 	for (i = num_state-1,inew = new_num-1; i >= 0; i--)
229 		if (forg[i] > ZERO_FREQ) {
230 			for (j = num_state-1, jnew = new_num-1; j >= 0; j--)
231 				if (forg[j] > ZERO_FREQ) {
232 					evec[i*num_state+j] = a[inew][jnew] / forg_sqrt[inew];
233 					inv_evec[i*num_state+j] = a[jnew][inew] * forg_sqrt[jnew];
234 					jnew--;
235 				} else {
236 					evec[i*num_state+j] = (i == j);
237 					inv_evec[i*num_state+j] = (i == j);
238 				}
239  			inew--;
240 		} else
241 		for (j=0; j < num_state; j++) {
242 			evec[i*num_state+j] = (i==j);
243 			inv_evec[i*num_state+j] = (i==j);
244 		}
245 
246     // Only print eigenvalues and eigenvectors if state space is manageable.
247 	if ((verbose_mode >= VB_MAX) && num_state < 30) {
248 		cout << "eigenvalues:";
249 		for (i = 0; i < num_state; i++)
250 			cout << "\t" << eval[i];
251 		cout << endl;
252 		cout << "eigenvectors:" << endl;
253 		for (i = 0; i < num_state; i++) {
254 			for (j = 0; j < num_state; j++) {
255 				cout << "\t" << evec[i*num_state+j];
256 			}
257 			cout << endl;
258 		}
259 		cout << "inv_eigenvectors:" << endl;
260 		for (i = 0; i < num_state; i++) {
261 			for (j = 0; j < num_state; j++) {
262 				cout << "\t" << inv_evec[i*num_state+j];
263 			}
264 			cout << endl;
265 		}
266 		cout << endl;
267 	}
268 
269 	/* check eigenvalue equation */
270 	error = 0.0;
271 	for (j = 0; j < num_state; j++) {
272 		for (i = 0; i < num_state; i++) {
273 			for (k = 0, zero = 0.0; k < num_state; k++) zero += b[i][k] * evec[k*num_state+j];
274             zero -= eval[j] * evec[i*num_state+j];
275             error = max(error, fabs(zero));
276 		}
277 	}
278 	if (error >= 0.1) {
279 		cout.precision(5);
280         cout.unsetf(ios::fixed);
281 		cout << "\nWARNING: Eigensystem doesn't satisfy eigenvalue equation! (gap=" << error << ")" << endl;
282 
283 		cout << " State frequencies (might be un-normalized): " << new_num << " states freq > " << ZERO_FREQ << endl;
284         double sum = 0.0;
285         cout.precision(7);
286 		for (i = 0; i < num_state; i++) {
287             cout << state_freq[i] << " ";
288             sum += state_freq[i];
289         }
290 		cout << endl;
291         cout << "sum = " << sum << endl;
292 		ASSERT(0);
293 	}
294 
295 	for (i=num_state-1; i>= 0; i--)
296 		delete [] b[i];
297 
298 	for (i=num_state-1; i>= 0; i--)
299 		delete [] a[i];
300 
301 	delete [] b;
302 	delete [] a;
303 	delete [] eval_new;
304 	delete [] off_diag;
305 	delete [] forg_sqrt;
306 	delete [] new_forg;
307 	delete [] forg;
308 
309 } // eigensystem_new
310 
311 
eigensystem_nonrev(double * rate_matrix,double * state_freq,double * eval,double * eval_imag,double * evec,double * inv_evec,int num_state)312 void EigenDecomposition::eigensystem_nonrev(
313 	double *rate_matrix, double *state_freq, double *eval, double *eval_imag,
314 	double *evec, double *inv_evec, int num_state)
315 {
316 	double *forg = new double[num_state];
317 	double *evali = new double[num_state];
318 	double *new_forg = new double[num_state];
319 	double *eval_new = new double[num_state];
320 	double **a = (double**)new double[num_state];
321 	double **b = (double**)new double[num_state];
322 	double **evec_new = (double**)new double[num_state];
323 	double **inv_evec_new = (double**)new double[num_state];
324 	int *ordr = new int[num_state + 1];
325 	int i, j, error, new_num, inew, jnew;
326 //	double zero;
327 
328 
329 	for (i=0; i < num_state; i++)
330 		a[i] = new double[num_state];
331 	for (i=0; i < num_state; i++)
332 		b[i] = new double[num_state];
333 	for (i=0; i < num_state; i++)
334 		evec_new[i] = new double[num_state];
335 	for (i=0; i < num_state; i++)
336 		inv_evec_new[i] = new double[num_state];
337 
338 	/* get relative transition matrix and frequencies */
339 	memcpy(forg, state_freq, num_state * sizeof(double));
340     // BQM 2015-09-07: normalize state frequencies to 1
341     double sum = 0.0;
342     for (i = 0; i < num_state; i++)
343         sum += forg[i];
344     sum = 1.0/sum;
345     for (i = 0; i < num_state; i++)
346         forg[i] *= sum;
347 
348 	for (i = 0; i < num_state; i++)
349 		memcpy(a[i], &rate_matrix[i*num_state], num_state * sizeof(double));
350 
351 	//rtfdata(a, forg, num_state);
352 	//    write (a, forg);
353 
354 	computeRateMatrix(a, forg, num_state); /* make 1 PAM rate matrix */
355 
356 	/* copy a to b */
357 	for (i = 0; i < num_state; i++)
358 		for (j = 0; j < num_state; j++)
359 			b[i][j] = a[i][j];
360 
361 	eliminateZero(b, forg, num_state, a, new_forg, new_num);
362 
363 	elmhes(a, ordr, new_num); /* compute eigenvalues and eigenvectors */
364 	//    writeInt (ordr);
365 
366 	eltran(a, evec_new, ordr, new_num);
367 
368 	//  writeMat (evec);
369 
370 	hqr2(new_num, 1, new_num, a, evec_new, eval_new, evali);
371 
372     // check that complex eigenvalues are conjugated
373     for (i = 0; i < new_num; i++)
374         if (evali[i] != 0.0) {
375             ASSERT(i < new_num-1 && evali[i+1] != 0.0);
376             i++;
377         }
378 
379 	luinverse(evec_new, inv_evec_new, num_state); /* inverse eigenvectors are in Ievc */
380 
381 	// now get back eigen
382 	//for (i = 0,inew = 0; i < num_state; i++)
383 	for (i = num_state-1,inew = new_num-1; i >= 0; i--)
384         if (forg[i] > ZERO_FREQ) {
385             eval[i] = eval_new[inew];
386             eval_imag[i] = evali[inew];
387             inew--;
388         } else {
389             eval[i] = 0.0;
390             eval_imag[i] = 0.0;
391         }
392 //		eval[i] = (forg[i] > ZERO) ? eval_new[inew--] : 0;
393 		//eval[i] = (forg[i] > ZERO) ? eval_new[inew++] : 0;
394 
395 	// calculate the actual eigenvectors of Q and its inverse matrix
396 	//for (i = 0, inew = 0; i < num_state; i++)
397 	for (i = num_state-1,inew = new_num-1; i >= 0; i--)
398 		if (forg[i] > ZERO_FREQ) {
399 // 			for (j = 0, jnew = 0; j < num_state; j++)
400 			for (j = num_state-1, jnew = new_num-1; j >= 0; j--)
401 				if (forg[j] > ZERO_FREQ) {
402 					evec[i*num_state+j] = evec_new[inew][jnew];
403 					inv_evec[i*num_state+j] = inv_evec_new[inew][jnew];
404 					//jnew++;
405 					jnew--;
406 				} else {
407 					evec[i*num_state+j] = (i == j);
408 					inv_evec[i*num_state+j] = (i == j);
409 				}
410 // 			inew++;
411  			inew--;
412 		} else
413 		for (j=0; j < num_state; j++) {
414 			evec[i*num_state+j] = (i==j);
415 			inv_evec[i*num_state+j] = (i==j);
416 		}
417 
418 	/* check eigenvalue equation */
419 	error = 0;
420 //	for (j = 0; j < num_state; j++) {
421 //		for (i = 0, zero = 0.0; i < num_state; i++) {
422 //			for (k = 0; k < num_state; k++) zero += b[i][k] * evec[k][j];
423 //			zero -= eval[j] * evec[i][j];
424 //			if (fabs(zero) > 1.0e-5) {
425 //				error = 1;
426 //				break;
427 //			}
428 //		}
429 //	}
430 	if (error) {
431 		cout << "\nWARNING: Eigensystem doesn't satisfy eigenvalue equation!\n";
432 		cout << "Rate matrix R: " << endl;
433 		for (i = 0; i < num_state; i++) {
434 			for (j = 0; j < num_state; j++) cout << rate_matrix[i*num_state+j] << " ";
435 			cout << endl;
436 		}
437 		cout << "State frequencies: " << endl;
438 		for (i = 0; i < num_state; i++) cout << state_freq[i] << " ";
439 		cout << endl;
440 	}
441 
442 	for (i=num_state-1; i>= 0; i--)
443 		delete [] inv_evec_new[i];
444 	for (i=num_state-1; i>= 0; i--)
445 		delete [] evec_new[i];
446 	for (i=num_state-1; i>= 0; i--)
447 		delete [] b[i];
448 	for (i=num_state-1; i>= 0; i--)
449 		delete [] a[i];
450 	delete [] ordr;
451 	delete [] inv_evec_new;
452 	delete [] evec_new;
453 	delete [] b;
454 	delete [] a;
455 	delete [] eval_new;
456 	delete [] new_forg;
457 	delete [] evali;
458 	delete [] forg;
459 
460 	//checkevector(evec, inv_evec, num_state); /* check whether inversion was OK */
461 
462 } /* eigensystem */
463 
464 
~EigenDecomposition()465 EigenDecomposition::~EigenDecomposition()
466 {
467 }
468 
469 /* make rate matrix with 0.01 expected substitutions per unit time */
computeRateMatrix(double ** a,double * stateFrqArr_,int num_state)470 void EigenDecomposition::computeRateMatrix(double **a, double *stateFrqArr_, int num_state) {
471 
472 /*
473 	if (myrate.isNsSyHeterogenous())
474 		return;
475 */
476 	int i, j;
477 	double delta, temp, sum;
478 	double *m = new double[num_state];
479 
480 	if (!ignore_state_freq)
481 	for (i = 0; i < num_state; i++) {
482 		for (j = 0; j < num_state; j++) {
483 			a[i][j] = stateFrqArr_[j]*a[i][j];
484 		}
485 	}
486 
487 	for (i = 0, sum = 0.0; i < num_state; i++) {
488 		for (j = 0, temp = 0.0; j < num_state; j++)
489             if (j != i)
490 			temp += a[i][j];
491 		m[i] = temp; /* row sum */
492 		sum += temp*stateFrqArr_[i]; /* exp. rate */
493 	}
494 
495 	if (normalize_matrix) {
496 		delta = total_num_subst / sum; /* 0.01 subst. per unit time */
497 
498 		for (i = 0; i < num_state; i++) {
499 			for (j = 0; j < num_state; j++) {
500 				if (i != j)
501 					a[i][j] = delta * a[i][j];
502 				else
503 					a[i][j] = delta * (-m[i]);
504 			}
505 		}
506 	} else {
507 		for (i = 0; i < num_state; i++)
508 			a[i][i] = -m[i];
509 	}
510 	delete [] m;
511 } /* onepamratematrix */
512 
eliminateZero(double ** mat,double * forg,int num,double ** new_mat,double * new_forg,int & new_num)513 void EigenDecomposition::eliminateZero(double **mat, double *forg, int num,
514 	double **new_mat, double *new_forg, int &new_num) {
515 	int i, j, inew, jnew;
516 	new_num = 0;
517 	for (i = 0; i < num; i++) {
518 		if (forg[i] > ZERO_FREQ)
519 			new_forg[new_num++] = forg[i];
520     }
521 	if (new_num == num) return;
522 	//writeDouble(forg, num);
523 	//writeMat(mat, num);
524 	for (i = 0, inew = 0; i < num; i++)
525 		if (forg[i] > ZERO_FREQ) {
526 			for (j = 0, jnew = 0; j < num; j++)
527 				if (forg[j] > ZERO_FREQ) {
528 					new_mat[inew][jnew] = mat[i][j];
529 					jnew++;
530 				}
531 			inew++;
532 		} else {
533             for (j = 0; j < num; j++)
534                 mat[i][j] = 0.0;
535             for (j = 0; j < num; j++)
536                 mat[j][i] = 0.0;
537         }
538 	if (verbose_mode >= VB_MED) {
539 		cout << "new_num_states = " << new_num << endl;
540         for (i = 0; i < new_num; i++) {
541             cout << new_mat[i][i] << " ";
542         }
543         cout << endl;
544     }
545 	//writeMat(new_mat, new_num);
546 	//writeDouble(new_forg, new_num);
547 }
548 
symmetrizeRateMatrix(double ** a,double * stateFrq,double * stateFrq_sqrt,int num_state)549 void EigenDecomposition::symmetrizeRateMatrix(double **a, double *stateFrq, double *stateFrq_sqrt, int num_state) {
550 	int i, j;
551 
552 	for (i = 0; i < num_state; i++)
553 		stateFrq_sqrt[i] = sqrt(stateFrq[i]);
554 	for (i = 0; i < num_state; i++) {
555         double tmp = 1.0/stateFrq_sqrt[i];
556 		for (j = 0; j < i; j++) {
557             a[j][i] *= stateFrq_sqrt[j]*tmp;
558             a[i][j] = a[j][i];
559 
560 //			a[i][j] *= stateFrq_sqrt[i] / stateFrq_sqrt[j];
561 //            a[j][i] = a[i][j];
562 
563 //            a[j][i] *= stateFrq_sqrt[j] / stateFrq_sqrt[i];
564 //			if (fabs(a[j][i] - a[i][j]) > 1e-5) {
565 //                cout << a[i][j] << "  " << a[j][i];
566 //                assert(0);
567 //            }
568 		}
569     }
570 }
571 
572 
tred2(double ** a,int n,double * d,double * e)573 void EigenDecomposition::tred2(double **a, int n, double *d, double *e)
574 {
575 	int l,k,j,i;
576 	double scale,hh,h,g,f;
577 
578 	for (i=n-1;i>0;i--) {
579 		l=i-1;
580 		h=scale=0.0;
581 		if (l > 0) {
582 			for (k=0;k<=l;k++)
583 				scale += fabs(a[i][k]);
584 			if (scale == 0.0)
585 				e[i]=a[i][l];
586 			else {
587 				for (k=0;k<=l;k++) {
588 					a[i][k] /= scale;
589 					h += a[i][k]*a[i][k];
590 				}
591 				f=a[i][l];
592 				g=(f >= 0.0 ? -sqrt(h) : sqrt(h));
593 				e[i]=scale*g;
594 				h -= f*g;
595 				a[i][l]=f-g;
596 				f=0.0;
597 				for (j=0;j<=l;j++) {
598 					a[j][i]=a[i][j]/h;
599 					g=0.0;
600 					for (k=0;k<=j;k++)
601 						g += a[j][k]*a[i][k];
602 					for (k=j+1;k<=l;k++)
603 						g += a[k][j]*a[i][k];
604 					e[j]=g/h;
605 					f += e[j]*a[i][j];
606 				}
607 				hh=f/(h+h);
608 				for (j=0;j<=l;j++) {
609 					f=a[i][j];
610 					e[j]=g=e[j]-hh*f;
611 					for (k=0;k<=j;k++)
612 						a[j][k] -= (f*e[k]+g*a[i][k]);
613 				}
614 			}
615 		} else
616 			e[i]=a[i][l];
617 		d[i]=h;
618 	}
619 	d[0]=0.0;
620 	e[0]=0.0;
621 	/* Contents of this loop can be omitted if eigenvectors not
622 			wanted except for statement d[i]=a[i][i]; */
623 	for (i=0;i<n;i++) {
624 		l=i;
625 		if (d[i] != 0.0) {
626 			for (j=0;j<l;j++) {
627 				g=0.0;
628 				for (k=0;k<l;k++)
629 					g += a[i][k]*a[k][j];
630 				for (k=0;k<l;k++)
631 					a[k][j] -= g*a[k][i];
632 			}
633 		}
634 		d[i]=a[i][i];
635 		a[i][i]=1.0;
636 		for (j=0;j<l;j++) a[j][i]=a[i][j]=0.0;
637 	}
638 }
639 
640 /**
641 	return a^2
642 */
sqr(double a)643 inline double sqr(double a) {
644 	return (a == 0.0) ? 0.0 : a*a;
645 }
646 
pythag(double a,double b)647 double pythag(double a, double b)
648 {
649 	double absa,absb;
650 	absa=fabs(a);
651 	absb=fabs(b);
652 	if (absa > absb) return absa*sqrt(1.0+sqr(absb/absa));
653 	else return (absb == 0.0 ? 0.0 : absb*sqrt(1.0+sqr(absa/absb)));
654 }
655 
656 #define NRANSI
657 #define SIGN(a,b) ((b) >= 0.0 ? fabs(a) : -fabs(a))
658 
tqli(double * d,double * e,int n,double ** z)659 void EigenDecomposition::tqli(double *d, double *e, int n, double **z)
660 {
661 	int m,l,iter,i,k;
662 	double s,r,p,g,f,dd,c,b;
663 
664 	for (i=1;i<n;i++) e[i-1]=e[i];
665 	e[n-1]=0.0;
666 	for (l=0;l<n;l++) {
667 		iter=0;
668 		do {
669 			for (m=l;m<n-1;m++) {
670 				dd=fabs(d[m])+fabs(d[m+1]);
671 				if ((double)(fabs(e[m])+dd) == dd) break;
672 			}
673 			if (m != l) {
674                 if (iter++ == 100) {
675 					outWarning("Too many iterations in tqli");
676                     break;
677                 }
678 
679 				g=(d[l+1]-d[l])/(2.0*e[l]);
680 				r=pythag(g,1.0);
681 				g=d[m]-d[l]+e[l]/(g+SIGN(r,g));
682 				s=c=1.0;
683 				p=0.0;
684 				for (i=m-1;i>=l;i--) {
685 					f=s*e[i];
686 					b=c*e[i];
687 					e[i+1]=(r=pythag(f,g));
688 					if (r == 0.0) {
689 						d[i+1] -= p;
690 						e[m]=0.0;
691 						break;
692 					}
693 					s=f/r;
694 					c=g/r;
695 					g=d[i+1]-p;
696 					r=(d[i]-g)*s+2.0*c*b;
697 					d[i+1]=g+(p=s*r);
698 					g=c*r-b;
699 					for (k=0;k<n;k++) {
700 						f=z[k][i+1];
701 						z[k][i+1]=s*z[k][i]+c*f;
702 						z[k][i]=c*z[k][i]-s*f;
703 					}
704 				}
705 				if (r == 0.0 && i >= l) continue;
706 				d[l] -= p;
707 				e[l]=g;
708 				e[m]=0.0;
709 			}
710 		} while (m != l);
711 	}
712 }
713 #undef SIGN
714 #undef NRANSI
715 
elmhes(double ** a,int * ordr,int n)716 void EigenDecomposition::elmhes(double **a, int *ordr, int n) {
717 	int m, j, i;
718 	double y, x;
719 
720 
721 	for (i = 0; i < n; i++)
722 		ordr[i] = 0;
723 	for (m = 2; m < n; m++) {
724 		x = 0.0;
725 		i = m;
726 		for (j = m; j <= n; j++) {
727 			if (fabs(a[j - 1][m - 2]) > fabs(x)) {
728 				x = a[j - 1][m - 2];
729 				i = j;
730 			}
731 		}
732 		ordr[m - 1] = i;      /* vector */
733 
734 		if (i != m) {
735 			for (j = m - 2; j < n; j++) {
736 				y = a[i - 1][j];
737 				a[i - 1][j] = a[m - 1][j];
738 				a[m - 1][j] = y;
739 			}
740 			for (j = 0; j < n; j++) {
741 				y = a[j][i - 1];
742 				a[j][i - 1] = a[j][m - 1];
743 				a[j][m - 1] = y;
744 			}
745 		}
746 		if (x != 0.0) {
747 			for (i = m; i < n; i++) {
748 				y = a[i][m - 2];
749 				if (y != 0.0) {
750 					y /= x;
751 					a[i][m - 2] = y;
752 					for (j = m - 1; j < n; j++)
753 						a[i][j] -= y * a[m - 1][j];
754 					for (j = 0; j < n; j++)
755 						a[j][m - 1] += y * a[j][i];
756 				}
757 			}
758 		}
759 	}
760 } /* elmhes */
761 
762 
eltran(double ** a,double ** zz,int * ordr,int n)763 void EigenDecomposition::eltran(double **a, double **zz, int *ordr, int n) {
764 	int i, j, m;
765 
766 
767 	for (i = 0; i < n; i++) {
768 		for (j = i + 1; j < n; j++) {
769 			zz[i][j] = 0.0;
770 			zz[j][i] = 0.0;
771 		}
772 		zz[i][i] = 1.0;
773 	}
774 	if (n <= 2)
775 		return;
776 	for (m = n - 1; m >= 2; m--) {
777 		for (i = m; i < n; i++)
778 			zz[i][m - 1] = a[i][m - 2];
779 		i = ordr[m - 1];
780 		if (i != m) {
781 			for (j = m - 1; j < n; j++) {
782 				zz[m - 1][j] = zz[i - 1][j];
783 				zz[i - 1][j] = 0.0;
784 			}
785 			zz[i - 1][m - 1] = 1.0;
786 		}
787 	}
788 } /* eltran */
789 
790 
mcdiv(double ar,double ai,double br,double bi,double * cr,double * ci)791 void EigenDecomposition::mcdiv(double ar, double ai, double br, double bi,
792                   double *cr, double *ci) {
793 	double s, ars, ais, brs, bis;
794 
795 
796 	s = fabs(br) + fabs(bi);
797 	ars = ar / s;
798 	ais = ai / s;
799 	brs = br / s;
800 	bis = bi / s;
801 	s = brs * brs + bis * bis;
802 	*cr = (ars * brs + ais * bis) / s;
803 	*ci = (ais * brs - ars * bis) / s;
804 } /* mcdiv */
805 
806 
hqr2(int n,int low,int hgh,double ** h,double ** zz,double * wr,double * wi)807 void EigenDecomposition::hqr2(int n, int low, int hgh, double **h,
808                  double **zz, double *wr, double *wi) {
809 	int i, j, k, l=0, m, en, na, itn, its;
810 	double p=0, q=0, r=0, s=0, t, w, x=0, y, ra, sa, vi, vr, z=0, norm, tst1, tst2;
811 	int notlas; /* boolean */
812 
813 
814 	norm = 0.0;
815 	k = 1;
816 	/* store isolated roots and compute matrix norm */
817 	for (i = 0; i < n; i++) {
818 		for (j = k - 1; j < n; j++)
819 			norm += fabs(h[i][j]);
820 		k = i + 1;
821 		if (i + 1 < low || i + 1 > hgh) {
822 			wr[i] = h[i][i];
823 			wi[i] = 0.0;
824 		}
825 	}
826 	en = hgh;
827 	t = 0.0;
828 	itn = n * 30;
829 	while (en >= low) {    /* search for next eigenvalues */
830 		its = 0;
831 		na = en - 1;
832 		while (en >= 1) {
833 			/* look for single small sub-diagonal element */
834 			for (l = en; l > low; l--) {
835 				s = fabs(h[l - 2][l - 2]) + fabs(h[l - 1][l - 1]);
836 
837 				if (s == 0.0)
838 					s = norm;
839 				tst1 = s;
840 				tst2 = tst1 + fabs(h[l - 1][l - 2]);
841 				if (tst2 == tst1)
842 					goto L100;
843 			}
844 			l = low;
845 		L100:
846 			x = h[en - 1][en - 1];    /* form shift */
847 			if (l == en || l == na)
848 				break;
849 			if (itn == 0) {
850 				/* all eigenvalues have not converged */
851 				cout << "\n\n\nHALT: PLEASE REPORT ERROR B TO DEVELOPERS\n\n\n";
852 				exit(1);
853 			}
854 			y = h[na - 1][na - 1];
855 			w = h[en - 1][na - 1] * h[na - 1][en - 1];
856 			/* form exceptional shift */
857 			if (its == 10 || its == 20) {
858 				t += x;
859 				for (i = low - 1; i < en; i++)
860 					h[i][i] -= x;
861 				s = fabs(h[en - 1][na - 1]) + fabs(h[na - 1][en - 3]);
862 				x = 0.75 * s;
863 				y = x;
864 				w = -0.4375 * s * s;
865 			}
866 			its++;
867 			itn--;
868 			/* look for two consecutive small sub-diagonal elements */
869 			for (m = en - 2; m >= l; m--) {
870 				z = h[m - 1][m - 1];
871 				r = x - z;
872 				s = y - z;
873 				p = (r * s - w) / h[m][m - 1] + h[m - 1][m];
874 				q = h[m][m] - z - r - s;
875 				r = h[m + 1][m];
876 				s = fabs(p) + fabs(q) + fabs(r);
877 				p /= s;
878 				q /= s;
879 				r /= s;
880 				if (m == l)
881 					break;
882 				tst1 = fabs(p) *
883 				       (fabs(h[m - 2][m - 2]) + fabs(z) + fabs(h[m][m]));
884 				tst2 = tst1 + fabs(h[m - 1][m - 2]) * (fabs(q) + fabs(r));
885 				if (tst2 == tst1)
886 					break;
887 			}
888 			for (i = m + 2; i <= en; i++) {
889 				h[i - 1][i - 3] = 0.0;
890 				if (i != m + 2)
891 					h[i - 1][i - 4] = 0.0;
892 			}
893 			for (k = m; k <= na; k++) {
894 				notlas = (k != na);
895 				if (k != m) {
896 					p = h[k - 1][k - 2];
897 					q = h[k][k - 2];
898 					r = 0.0;
899 					if (notlas)
900 						r = h[k + 1][k - 2];
901 					x = fabs(p) + fabs(q) + fabs(r);
902 					if (x != 0.0) {
903 						p /= x;
904 						q /= x;
905 						r /= x;
906 					}
907 				}
908 				if (x != 0.0) {
909 					if (p < 0.0) /* sign */
910 						s = - sqrt(p * p + q * q + r * r);
911 					else
912 						s = sqrt(p * p + q * q + r * r);
913 					if (k != m)
914 						h[k - 1][k - 2] = -s * x;
915 					else {
916 						if (l != m)
917 							h[k - 1][k - 2] = -h[k - 1][k - 2];
918 					}
919 					p += s;
920 					x = p / s;
921 					y = q / s;
922 					z = r / s;
923 					q /= p;
924 					r /= p;
925 					if (!notlas) {
926 						for (j = k - 1; j < n; j++) {    /* row modification */
927 							p = h[k - 1][j] + q * h[k][j];
928 							h[k - 1][j] -= p * x;
929 							h[k][j] -= p * y;
930 						}
931 						j = (en < (k + 3)) ? en : (k + 3); /* min */
932 						for (i = 0; i < j; i++) {    /* column modification */
933 							p = x * h[i][k - 1] + y * h[i][k];
934 							h[i][k - 1] -= p;
935 							h[i][k] -= p * q;
936 						}
937 						/* accumulate transformations */
938 						for (i = low - 1; i < hgh; i++) {
939 							p = x * zz[i][k - 1] + y * zz[i][k];
940 							zz[i][k - 1] -= p;
941 							zz[i][k] -= p * q;
942 						}
943 					} else {
944 						for (j = k - 1; j < n; j++) {    /* row modification */
945 							p = h[k - 1][j] + q * h[k][j] + r * h[k + 1][j];
946 							h[k - 1][j] -= p * x;
947 							h[k][j] -= p * y;
948 							h[k + 1][j] -= p * z;
949 						}
950 						j = (en < (k + 3)) ? en : (k + 3); /* min */
951 						for (i = 0; i < j; i++) {    /* column modification */
952 							p = x * h[i][k - 1] + y * h[i][k] + z * h[i][k + 1];
953 							h[i][k - 1] -= p;
954 							h[i][k] -= p * q;
955 							h[i][k + 1] -= p * r;
956 						}
957 						/* accumulate transformations */
958 						for (i = low - 1; i < hgh; i++) {
959 							p = x * zz[i][k - 1] + y * zz[i][k] +
960 							    z * zz[i][k + 1];
961 							zz[i][k - 1] -= p;
962 							zz[i][k] -= p * q;
963 							zz[i][k + 1] -= p * r;
964 						}
965 					}
966 				}
967 			}           /* for k */
968 		}               /* while infinite loop */
969 		if (l == en) {           /* one root found */
970 			h[en - 1][en - 1] = x + t;
971 			wr[en - 1] = h[en - 1][en - 1];
972 			wi[en - 1] = 0.0;
973 			en = na;
974 			continue;
975 		}
976 		y = h[na - 1][na - 1];
977 		w = h[en - 1][na - 1] * h[na - 1][en - 1];
978 		p = (y - x) / 2.0;
979 		q = p * p + w;
980 		z = sqrt(fabs(q));
981 		h[en - 1][en - 1] = x + t;
982 		x = h[en - 1][en - 1];
983 		h[na - 1][na - 1] = y + t;
984 		if (q >= 0.0) {           /* real pair */
985 			if (p < 0.0) /* sign */
986 				z = p - fabs(z);
987 			else
988 				z = p + fabs(z);
989 			wr[na - 1] = x + z;
990 			wr[en - 1] = wr[na - 1];
991 			if (z != 0.0)
992 				wr[en - 1] = x - w / z;
993 			wi[na - 1] = 0.0;
994 			wi[en - 1] = 0.0;
995 			x = h[en - 1][na - 1];
996 			s = fabs(x) + fabs(z);
997 			p = x / s;
998 			q = z / s;
999 			r = sqrt(p * p + q * q);
1000 			p /= r;
1001 			q /= r;
1002 			for (j = na - 1; j < n; j++) {    /* row modification */
1003 				z = h[na - 1][j];
1004 				h[na - 1][j] = q * z + p * h[en - 1][j];
1005 				h[en - 1][j] = q * h[en - 1][j] - p * z;
1006 			}
1007 			for (i = 0; i < en; i++) {    /* column modification */
1008 				z = h[i][na - 1];
1009 				h[i][na - 1] = q * z + p * h[i][en - 1];
1010 				h[i][en - 1] = q * h[i][en - 1] - p * z;
1011 			}
1012 			/* accumulate transformations */
1013 			for (i = low - 1; i < hgh; i++) {
1014 				z = zz[i][na - 1];
1015 				zz[i][na - 1] = q * z + p * zz[i][en - 1];
1016 				zz[i][en - 1] = q * zz[i][en - 1] - p * z;
1017 			}
1018 		} else {           /* complex pair */
1019 			wr[na - 1] = x + p;
1020 			wr[en - 1] = x + p;
1021 			wi[na - 1] = z;
1022 			wi[en - 1] = -z;
1023 		}
1024 		en -= 2;
1025 	}                   /* while en >= low */
1026 	/* backsubstitute to find vectors of upper triangular form */
1027 	if (norm != 0.0) {
1028 		for (en = n; en >= 1; en--) {
1029 			p = wr[en - 1];
1030 			q = wi[en - 1];
1031 			na = en - 1;
1032 			if (q == 0.0) {/* real vector */
1033 				m = en;
1034 				h[en - 1][en - 1] = 1.0;
1035 				if (na != 0) {
1036 					for (i = en - 2; i >= 0; i--) {
1037 						w = h[i][i] - p;
1038 						r = 0.0;
1039 						for (j = m - 1; j < en; j++)
1040 							r += h[i][j] * h[j][en - 1];
1041 						if (wi[i] < 0.0) {
1042 							z = w;
1043 							s = r;
1044 						} else {
1045 							m = i + 1;
1046 							if (wi[i] == 0.0) {
1047 								t = w;
1048 								if (t == 0.0) {
1049 									tst1 = norm;
1050 									t = tst1;
1051 									do {
1052 										t = 0.01 * t;
1053 										tst2 = norm + t;
1054 									} while (tst2 > tst1);
1055 								}
1056 								h[i][en - 1] = -(r / t);
1057 							} else {    /* solve real equations */
1058 								x = h[i][i + 1];
1059 								y = h[i + 1][i];
1060 								q = (wr[i] - p) * (wr[i] - p) + wi[i] * wi[i];
1061 								t = (x * s - z * r) / q;
1062 								h[i][en - 1] = t;
1063 								if (fabs(x) > fabs(z))
1064 									h[i + 1][en - 1] = (-r - w * t) / x;
1065 								else
1066 									h[i + 1][en - 1] = (-s - y * t) / z;
1067 							}
1068 							/* overflow control */
1069 							t = fabs(h[i][en - 1]);
1070 							if (t != 0.0) {
1071 								tst1 = t;
1072 								tst2 = tst1 + 1.0 / tst1;
1073 								if (tst2 <= tst1) {
1074 									for (j = i; j < en; j++)
1075 										h[j][en - 1] /= t;
1076 								}
1077 							}
1078 						}
1079 					}
1080 				}
1081 			} else if (q > 0.0) {
1082 				m = na;
1083 				if (fabs(h[en - 1][na - 1]) > fabs(h[na - 1][en - 1])) {
1084 					h[na - 1][na - 1] = q / h[en - 1][na - 1];
1085 					h[na - 1][en - 1] = (p - h[en - 1][en - 1]) /
1086 					                    h[en - 1][na - 1];
1087 				} else
1088 					mcdiv(0.0, -h[na - 1][en - 1], h[na - 1][na - 1] - p, q,
1089 					      &h[na - 1][na - 1], &h[na - 1][en - 1]);
1090 				h[en - 1][na - 1] = 0.0;
1091 				h[en - 1][en - 1] = 1.0;
1092 				if (en != 2) {
1093 					for (i = en - 3; i >= 0; i--) {
1094 						w = h[i][i] - p;
1095 						ra = 0.0;
1096 						sa = 0.0;
1097 						for (j = m - 1; j < en; j++) {
1098 							ra += h[i][j] * h[j][na - 1];
1099 							sa += h[i][j] * h[j][en - 1];
1100 						}
1101 						if (wi[i] < 0.0) {
1102 							z = w;
1103 							r = ra;
1104 							s = sa;
1105 						} else {
1106 							m = i + 1;
1107 							if (wi[i] == 0.0)
1108 								mcdiv(-ra, -sa, w, q, &h[i][na - 1],
1109 								      &h[i][en - 1]);
1110 							else {    /* solve complex equations */
1111 								x = h[i][i + 1];
1112 								y = h[i + 1][i];
1113 								vr = (wr[i] - p) * (wr[i] - p);
1114 								vr = vr + wi[i] * wi[i] - q * q;
1115 								vi = (wr[i] - p) * 2.0 * q;
1116 								if (vr == 0.0 && vi == 0.0) {
1117 									tst1 = norm * (fabs(w) + fabs(q) + fabs(x) +
1118 									               fabs(y) + fabs(z));
1119 									vr = tst1;
1120 									do {
1121 										vr = 0.01 * vr;
1122 										tst2 = tst1 + vr;
1123 									} while (tst2 > tst1);
1124 								}
1125 								mcdiv(x * r - z * ra + q * sa,
1126 								      x * s - z * sa - q * ra, vr, vi,
1127 								      &h[i][na - 1], &h[i][en - 1]);
1128 								if (fabs(x) > fabs(z) + fabs(q)) {
1129 									h[i + 1]
1130 									[na - 1] = (q * h[i][en - 1] -
1131 									            w * h[i][na - 1] - ra) / x;
1132 									h[i + 1][en - 1] = (-sa - w * h[i][en - 1] -
1133 									                    q * h[i][na - 1]) / x;
1134 								} else
1135 									mcdiv(-r - y * h[i][na - 1],
1136 									      -s - y * h[i][en - 1], z, q,
1137 									      &h[i + 1][na - 1], &h[i + 1][en - 1]);
1138 							}
1139 							/* overflow control */
1140 							t = (fabs(h[i][na - 1]) > fabs(h[i][en - 1])) ?
1141 							    fabs(h[i][na - 1]) : fabs(h[i][en - 1]);
1142 							if (t != 0.0) {
1143 								tst1 = t;
1144 								tst2 = tst1 + 1.0 / tst1;
1145 								if (tst2 <= tst1) {
1146 									for (j = i; j < en; j++) {
1147 										h[j][na - 1] /= t;
1148 										h[j][en - 1] /= t;
1149 									}
1150 								}
1151 							}
1152 						}
1153 					}
1154 				}
1155 			}
1156 		}
1157 		/* end back substitution. vectors of isolated roots */
1158 		for (i = 0; i < n; i++) {
1159 			if (i + 1 < low || i + 1 > hgh) {
1160 				for (j = i; j < n; j++)
1161 					zz[i][j] = h[i][j];
1162 			}
1163 		}
1164 		/* multiply by transformation matrix to give vectors of
1165 		 * original full matrix. */
1166 		for (j = n - 1; j >= low - 1; j--) {
1167 			m = ((j + 1) < hgh) ? (j + 1) : hgh; /* min */
1168 			for (i = low - 1; i < hgh; i++) {
1169 				z = 0.0;
1170 				for (k = low - 1; k < m; k++)
1171 					z += zz[i][k] * h[k][j];
1172 				zz[i][j] = z;
1173 			}
1174 		}
1175 	}
1176 	return;
1177 } /* hqr2 */
1178 
luinverse(double ** inmat,double ** imtrx,int size)1179 void EigenDecomposition::luinverse(double **inmat, double **imtrx, int size) {
1180 	double eps = 1.0e-20; /* ! */
1181 	int i, j, k, l, maxi=0, idx, ix, jx;
1182 	double sum, tmp, maxb, aw;
1183 	int *index = new int[size];
1184 	double *wk;
1185 	double **omtrx = (double**) new double[size];
1186 
1187 	for (i = 0; i < size; i++)
1188 		omtrx[i] = new double[size];
1189 
1190 	/* copy inmat to omtrx */
1191 	for (i = 0; i < size; i++)
1192 		for (j = 0; j < size; j++)
1193 			omtrx[i][j] = inmat[i][j];
1194 
1195 	wk = (double *) calloc((size_t)size, sizeof(double));
1196 	aw = 1.0;
1197 	for (i = 0; i < size; i++) {
1198 		maxb = 0.0;
1199 		for (j = 0; j < size; j++) {
1200 			if (fabs(omtrx[i][j]) > maxb)
1201 				maxb = fabs(omtrx[i][j]);
1202 		}
1203 		if (maxb == 0.0) {
1204 			/* Singular matrix */
1205 			ASSERT(0 && "\n\n\nHALT: PLEASE REPORT ERROR C TO DEVELOPERS\n\n\n");
1206 		}
1207 		wk[i] = 1.0 / maxb;
1208 	}
1209 	for (j = 0; j < size; j++) {
1210 		for (i = 0; i < j; i++) {
1211 			sum = omtrx[i][j];
1212 			for (k = 0; k < i; k++)
1213 				sum -= omtrx[i][k] * omtrx[k][j];
1214 			omtrx[i][j] = sum;
1215 		}
1216 		maxb = 0.0;
1217 		for (i = j; i < size; i++) {
1218 			sum = omtrx[i][j];
1219 			for (k = 0; k < j; k++)
1220 				sum -= omtrx[i][k] * omtrx[k][j];
1221 			omtrx[i][j] = sum;
1222 			tmp = wk[i] * fabs(sum);
1223 			if (tmp >= maxb) {
1224 				maxb = tmp;
1225 				maxi = i;
1226 			}
1227 		}
1228 		if (j != maxi) {
1229 			for (k = 0; k < size; k++) {
1230 				tmp = omtrx[maxi][k];
1231 				omtrx[maxi][k] = omtrx[j][k];
1232 				omtrx[j][k] = tmp;
1233 			}
1234 			aw = -aw;
1235 			wk[maxi] = wk[j];
1236 		}
1237 		index[j] = maxi;
1238 		if (omtrx[j][j] == 0.0)
1239 			omtrx[j][j] = eps;
1240 		if (j != size - 1) {
1241 			tmp = 1.0 / omtrx[j][j];
1242 			for (i = j + 1; i < size; i++)
1243 				omtrx[i][j] *= tmp;
1244 		}
1245 	}
1246 	for (jx = 0; jx < size; jx++) {
1247 		for (ix = 0; ix < size; ix++)
1248 			wk[ix] = 0.0;
1249 		wk[jx] = 1.0;
1250 		l = -1;
1251 		for (i = 0; i < size; i++) {
1252 			idx = index[i];
1253 			sum = wk[idx];
1254 			wk[idx] = wk[i];
1255 			if (l != -1) {
1256 				for (j = l; j < i; j++)
1257 					sum -= omtrx[i][j] * wk[j];
1258 			} else if (sum != 0.0)
1259 				l = i;
1260 			wk[i] = sum;
1261 		}
1262 		for (i = size - 1; i >= 0; i--) {
1263 			sum = wk[i];
1264 			for (j = i + 1; j < size; j++)
1265 				sum -= omtrx[i][j] * wk[j];
1266 			wk[i] = sum / omtrx[i][i];
1267 		}
1268 		for (ix = 0; ix < size; ix++)
1269 			imtrx[ix][jx] = wk[ix];
1270 	}
1271 	free((char *)wk);
1272 	wk = NULL;
1273 	for (i = size-1; i >= 0; i--)
1274 		delete [] omtrx[i];
1275 	delete [] omtrx;
1276 	delete [] index;
1277 } /* luinverse */
1278 
checkevector(double * evec,double * ivec,int nn)1279 void EigenDecomposition::checkevector(double *evec, double *ivec, int nn) {
1280 	int i, j, ia, ib, ic, error;
1281 	double **matx = (double**) new double [nn];
1282 	double sum;
1283 
1284 	for (i = 0; i < nn; i++)
1285 		matx[i] = new double[nn];
1286 
1287 	/* multiply matrix of eigenvectors and its inverse */
1288 	for (ia = 0; ia < nn; ia++) {
1289 		for (ic = 0; ic < nn; ic++) {
1290 			sum = 0.0;
1291 			for (ib = 0; ib < nn; ib++) sum += evec[ia*nn+ib] * ivec[ib*nn+ic];
1292 			matx[ia][ic] = sum;
1293 		}
1294 	}
1295 	/* check whether the unitary matrix is obtained */
1296 	error = 0;
1297 	for (i = 0; i < nn; i++) {
1298 		for (j = 0; j < nn; j++) {
1299 			if (i == j) {
1300 				if (fabs(matx[i][j] - 1.0) > 1.0e-5)
1301 					error = 1;
1302 			} else {
1303 				if (fabs(matx[i][j]) > 1.0e-5)
1304 					error = 1;
1305 			}
1306 		}
1307 	}
1308 	if (error) {
1309 		cout << "\nWARNING: Inversion of eigenvector matrix not perfect!\n";
1310 	}
1311 
1312 	for (i = nn-1; i >= 0; i--)
1313 		delete [] matx[i];
1314 	delete [] matx;
1315 } /* checkevector */
1316