Algorithms_in_C  1.0.0
Set of algorithms implemented in C.
qr_eigen_values.c File Reference

Compute real eigen values and eigen vectors of a symmetric matrix using QR decomposition method. More...

#include <assert.h>
#include <math.h>
#include <stdio.h>
#include <stdlib.h>
#include <time.h>
#include "qr_decompose.h"
Include dependency graph for qr_eigen_values.c:

Macros

#define LIMS   9
 
#define EPSILON   1e-10
 

Functions

void create_matrix (double **A, int N)
 
double ** mat_mul (double **A, double **B, double **OUT, int R1, int C1, int R2, int C2)
 
double eigen_values (double **A, double *eigen_vals, int mat_size, char debug_print)
 
void test1 ()
 
void test2 ()
 
int main (int argc, char **argv)
 

Detailed Description

Compute real eigen values and eigen vectors of a symmetric matrix using QR decomposition method.

Author
Krishna Vedala

Macro Definition Documentation

◆ EPSILON

#define EPSILON   1e-10

accuracy tolerance limit

◆ LIMS

#define LIMS   9

limit of range of matrix values

Function Documentation

◆ create_matrix()

void create_matrix ( double **  A,
int  N 
)

create a square matrix of given size with random elements

Parameters
[out]Amatrix to create (must be pre-allocated in memory)
[in]Nmatrix size
27 {
28  int i, j, tmp, lim2 = LIMS >> 1;
29 
30 #ifdef _OPENMP
31 #pragma omp for
32 #endif
33  for (i = 0; i < N; i++)
34  {
35  A[i][i] = (rand() % LIMS) - lim2;
36  for (j = i + 1; j < N; j++)
37  {
38  tmp = (rand() % LIMS) - lim2;
39  A[i][j] = tmp;
40  A[j][i] = tmp;
41  }
42  }
43 }

◆ eigen_values()

double eigen_values ( double **  A,
double *  eigen_vals,
int  mat_size,
char  debug_print 
)

Compute eigen values using iterative shifted QR decomposition algorithm as follows:

  1. Use last diagonal element of A as eigen value approximation \(c\)
  2. Shift diagonals of matrix \(A' = A - cI\)
  3. Decompose matrix \(A'=QR\)
  4. Compute next approximation \(A'_1 = RQ \)
  5. Shift diagonals back \(A_1 = A'_1 + cI\)
  6. Termination condition check: last element below diagonal is almost 0
    1. If not 0, go back to step 1 with the new approximation \(A_1\)
    2. If 0, continue to step 7
  7. Save last known \(c\) as the eigen value.
  8. Are all eigen values found?
    1. If not, remove last row and column of \(A_1\) and go back to step 1.
    2. If yes, stop.
Note
The matrix \(A\) gets modified
Parameters
[in,out]Amatrix to compute eigen values for
[out]eig_valsresultant vector containing computed eigen values
[in]mat_sizematrix size
[in]debug_print1 to print intermediate Q & R matrices, 0 for not to
Returns
time for computation in seconds
105 {
106  double **R = (double **)malloc(sizeof(double *) * mat_size);
107  double **Q = (double **)malloc(sizeof(double *) * mat_size);
108 
109  if (!eigen_vals)
110  {
111  perror("Output eigen value vector cannot be NULL!");
112  return -1;
113  }
114  else if (!Q || !R)
115  {
116  perror("Unable to allocate memory for Q & R!");
117  return -1;
118  }
119 
120  /* allocate dynamic memory for matrices */
121  for (int i = 0; i < mat_size; i++)
122  {
123  R[i] = (double *)malloc(sizeof(double) * mat_size);
124  Q[i] = (double *)malloc(sizeof(double) * mat_size);
125  if (!Q[i] || !R[i])
126  {
127  perror("Unable to allocate memory for Q & R.");
128  return -1;
129  }
130  }
131 
132  if (debug_print)
133  print_matrix(A, mat_size, mat_size);
134 
135  int rows = mat_size, columns = mat_size;
136  int counter = 0, num_eigs = rows - 1;
137  double last_eig = 0;
138 
139  clock_t t1 = clock();
140  while (num_eigs > 0) /* continue till all eigen values are found */
141  {
142  /* iterate with QR decomposition */
143  while (fabs(A[num_eigs][num_eigs - 1]) > EPSILON)
144  {
145  last_eig = A[num_eigs][num_eigs];
146  for (int i = 0; i < rows; i++) A[i][i] -= last_eig; /* A - cI */
147  qr_decompose(A, Q, R, rows, columns);
148 
149  if (debug_print)
150  {
151  print_matrix(A, rows, columns);
152  print_matrix(Q, rows, columns);
153  print_matrix(R, columns, columns);
154  printf("-------------------- %d ---------------------\n",
155  ++counter);
156  }
157 
158  mat_mul(R, Q, A, columns, columns, rows, columns);
159  for (int i = 0; i < rows; i++) A[i][i] += last_eig; /* A + cI */
160  }
161 
162  /* store the converged eigen value */
163  eigen_vals[num_eigs] = last_eig;
164 
165  if (debug_print)
166  {
167  printf("========================\n");
168  printf("Eigen value: % g,\n", last_eig);
169  printf("========================\n");
170  }
171 
172  num_eigs--;
173  rows--;
174  columns--;
175  }
176  eigen_vals[0] = A[0][0];
177  double dtime = (double)(clock() - t1) / CLOCKS_PER_SEC;
178 
179  if (debug_print)
180  {
181  print_matrix(R, mat_size, mat_size);
182  print_matrix(Q, mat_size, mat_size);
183  }
184 
185  /* cleanup dynamic memory */
186  for (int i = 0; i < mat_size; i++)
187  {
188  free(R[i]);
189  free(Q[i]);
190  }
191  free(R);
192  free(Q);
193 
194  return dtime;
195 }
Here is the call graph for this function:

◆ main()

int main ( int  argc,
char **  argv 
)

main function

297 {
298  srand(time(NULL));
299 
300  int mat_size = 5;
301  if (argc == 2)
302  mat_size = atoi(argv[1]);
303  else
304  { // if invalid input argument is given run tests
305  test1();
306  test2();
307  printf("Usage: ./qr_eigen_values [mat_size]\n");
308  return 0;
309  }
310 
311  if (mat_size < 2)
312  {
313  fprintf(stderr, "Matrix size should be > 2\n");
314  return -1;
315  }
316 
317  int i, rows = mat_size, columns = mat_size;
318 
319  double **A = (double **)malloc(sizeof(double *) * mat_size);
320  /* number of eigen values = matrix size */
321  double *eigen_vals = (double *)malloc(sizeof(double) * mat_size);
322  if (!eigen_vals)
323  {
324  perror("Unable to allocate memory for eigen values!");
325  return -1;
326  }
327  for (i = 0; i < mat_size; i++)
328  {
329  A[i] = (double *)malloc(sizeof(double) * mat_size);
330  }
331 
332  /* create a random matrix */
333  create_matrix(A, mat_size);
334 
335  print_matrix(A, mat_size, mat_size);
336 
337  double dtime = eigen_values(A, eigen_vals, mat_size, 0);
338  printf("Eigen vals: ");
339  for (i = 0; i < mat_size; i++) printf("% 9.4g\t", eigen_vals[i]);
340  printf("\nTime taken to compute: % .4g sec\n", dtime);
341 
342  for (int i = 0; i < mat_size; i++) free(A[i]);
343  free(A);
344  free(eigen_vals);
345  return 0;
346 }
Here is the call graph for this function:

◆ mat_mul()

double** mat_mul ( double **  A,
double **  B,
double **  OUT,
int  R1,
int  C1,
int  R2,
int  C2 
)

Perform multiplication of two matrices.

  • R2 must be equal to C1
  • Resultant matrix size should be R1xC2
    Parameters
    [in]Afirst matrix to multiply
    [in]Bsecond matrix to multiply
    [out]OUToutput matrix (must be pre-allocated)
    [in]R1number of rows of first matrix
    [in]C1number of columns of first matrix
    [in]R2number of rows of second matrix
    [in]C2number of columns of second matrix
    Returns
    pointer to resultant matrix
60 {
61  if (C1 != R2)
62  {
63  perror("Matrix dimensions mismatch!");
64  return OUT;
65  }
66 
67  int i;
68 #ifdef _OPENMP
69 #pragma omp for
70 #endif
71  for (i = 0; i < R1; i++)
72  for (int j = 0; j < C2; j++)
73  {
74  OUT[i][j] = 0.f;
75  for (int k = 0; k < C1; k++) OUT[i][j] += A[i][k] * B[k][j];
76  }
77  return OUT;
78 }

◆ test1()

void test1 ( )

test function to compute eigen values of a 2x2 matrix

\[\begin{bmatrix} 5 & 7\\ 7 & 11 \end{bmatrix}\]

which are approximately, {15.56158, 0.384227}

206 {
207  int mat_size = 2;
208  double X[][2] = {{5, 7}, {7, 11}};
209  double y[] = {15.56158, 0.384227}; // corresponding y-values
210  double eig_vals[2];
211 
212  // The following steps are to convert a "double[][]" to "double **"
213  double **A = (double **)malloc(mat_size * sizeof(double *));
214  for (int i = 0; i < mat_size; i++) A[i] = X[i];
215 
216  printf("------- Test 1 -------\n");
217 
218  double dtime = eigen_values(A, eig_vals, mat_size, 0);
219 
220  for (int i = 0; i < mat_size; i++)
221  {
222  printf("%d/5 Checking for %.3g --> ", i + 1, y[i]);
223  char result = 0;
224  for (int j = 0; j < mat_size && !result; j++)
225  {
226  if (fabs(y[i] - eig_vals[j]) < 0.1)
227  {
228  result = 1;
229  printf("(%.3g) ", eig_vals[j]);
230  }
231  }
232 
233  // ensure that i^th expected eigen value was computed
234  assert(result != 0);
235  printf("found\n");
236  }
237  printf("Test 1 Passed in %.3g sec\n\n", dtime);
238  free(A);
239 }
Here is the call graph for this function:

◆ test2()

void test2 ( )

test function to compute eigen values of a 2x2 matrix

\[\begin{bmatrix} -4& 4& 2& 0& -3\\ 4& -4& 4& -3& -1\\ 2& 4& 4& 3& -3\\ 0& -3& 3& -1&-1\\ -3& -1& -3& -3& 0 \end{bmatrix}\]

which are approximately, {9.27648, -9.26948, 2.0181, -1.03516, -5.98994}

253 {
254  int mat_size = 5;
255  double X[][5] = {{-4, 4, 2, 0, -3},
256  {4, -4, 4, -3, -1},
257  {2, 4, 4, 3, -3},
258  {0, -3, 3, -1, -3},
259  {-3, -1, -3, -3, 0}};
260  double y[] = {9.27648, -9.26948, 2.0181, -1.03516,
261  -5.98994}; // corresponding y-values
262  double eig_vals[5];
263 
264  // The following steps are to convert a "double[][]" to "double **"
265  double **A = (double **)malloc(mat_size * sizeof(double *));
266  for (int i = 0; i < mat_size; i++) A[i] = X[i];
267 
268  printf("------- Test 2 -------\n");
269 
270  double dtime = eigen_values(A, eig_vals, mat_size, 0);
271 
272  for (int i = 0; i < mat_size; i++)
273  {
274  printf("%d/5 Checking for %.3g --> ", i + 1, y[i]);
275  char result = 0;
276  for (int j = 0; j < mat_size && !result; j++)
277  {
278  if (fabs(y[i] - eig_vals[j]) < 0.1)
279  {
280  result = 1;
281  printf("(%.3g) ", eig_vals[j]);
282  }
283  }
284 
285  // ensure that i^th expected eigen value was computed
286  assert(result != 0);
287  printf("found\n");
288  }
289  printf("Test 2 Passed in %.3g sec\n\n", dtime);
290  free(A);
291 }
Here is the call graph for this function:
EPSILON
#define EPSILON
Definition: qr_eigen_values.c:19
test2
void test2()
Definition: qr_eigen_values.c:252
N
#define N
Definition: sol1.c:109
qr_decompose
void qr_decompose(double **A, double **Q, double **R, int M, int N)
Definition: qr_decompose.h:142
test1
void test1()
Definition: qr_eigen_values.c:205
print_matrix
void print_matrix(double **A, int M, int N)
Definition: qr_decompose.h:22
mat_mul
double ** mat_mul(double **A, double **B, double **OUT, int R1, int C1, int R2, int C2)
Definition: qr_eigen_values.c:58
LIMS
#define LIMS
Definition: qr_eigen_values.c:18
eigen_values
double eigen_values(double **A, double *eigen_vals, int mat_size, char debug_print)
Definition: qr_eigen_values.c:103
create_matrix
void create_matrix(double **A, int N)
Definition: qr_eigen_values.c:26