File: igraph_sparsemat5.c

package info (click to toggle)
igraph 0.7.1-2
  • links: PTS, VCS
  • area: main
  • in suites: jessie, jessie-kfreebsd
  • size: 19,180 kB
  • ctags: 16,601
  • sloc: ansic: 188,037; sh: 26,731; cpp: 18,275; yacc: 1,164; makefile: 959; lex: 484; xml: 378
file content (351 lines) | stat: -rw-r--r-- 10,826 bytes parent folder | download | duplicates (3)
1
2
3
4
5
6
7
8
9
10
11
12
13
14
15
16
17
18
19
20
21
22
23
24
25
26
27
28
29
30
31
32
33
34
35
36
37
38
39
40
41
42
43
44
45
46
47
48
49
50
51
52
53
54
55
56
57
58
59
60
61
62
63
64
65
66
67
68
69
70
71
72
73
74
75
76
77
78
79
80
81
82
83
84
85
86
87
88
89
90
91
92
93
94
95
96
97
98
99
100
101
102
103
104
105
106
107
108
109
110
111
112
113
114
115
116
117
118
119
120
121
122
123
124
125
126
127
128
129
130
131
132
133
134
135
136
137
138
139
140
141
142
143
144
145
146
147
148
149
150
151
152
153
154
155
156
157
158
159
160
161
162
163
164
165
166
167
168
169
170
171
172
173
174
175
176
177
178
179
180
181
182
183
184
185
186
187
188
189
190
191
192
193
194
195
196
197
198
199
200
201
202
203
204
205
206
207
208
209
210
211
212
213
214
215
216
217
218
219
220
221
222
223
224
225
226
227
228
229
230
231
232
233
234
235
236
237
238
239
240
241
242
243
244
245
246
247
248
249
250
251
252
253
254
255
256
257
258
259
260
261
262
263
264
265
266
267
268
269
270
271
272
273
274
275
276
277
278
279
280
281
282
283
284
285
286
287
288
289
290
291
292
293
294
295
296
297
298
299
300
301
302
303
304
305
306
307
308
309
310
311
312
313
314
315
316
317
318
319
320
321
322
323
324
325
326
327
328
329
330
331
332
333
334
335
336
337
338
339
340
341
342
343
344
345
346
347
348
349
350
351
/* -*- mode: C -*-  */
/* 
   IGraph library.
   Copyright (C) 2009-2012  Gabor Csardi <csardi.gabor@gmail.com>
   334 Harvard st, Cambridge MA, 02139 USA
   
   This program is free software; you can redistribute it and/or modify
   it under the terms of the GNU General Public License as published by
   the Free Software Foundation; either version 2 of the License, or
   (at your option) any later version.
   
   This program is distributed in the hope that it will be useful,
   but WITHOUT ANY WARRANTY; without even the implied warranty of
   MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the
   GNU General Public License for more details.
   
   You should have received a copy of the GNU General Public License
   along with this program; if not, write to the Free Software
   Foundation, Inc., 51 Franklin Street, Fifth Floor, Boston, MA 
   02110-1301 USA

*/

#include <igraph.h>
#include <igraph_sparsemat.h>

#define EPS 1e-13


/* Generic test for 1x1 matrices */
void test_1x1(igraph_real_t value) {
  igraph_sparsemat_t A, B;
  igraph_matrix_t values, vectors;
  igraph_vector_t values2;
  igraph_arpack_options_t options;

  igraph_arpack_options_init(&options);

  igraph_sparsemat_init(&A, 1, 1, 1);
  igraph_sparsemat_entry(&A, 0, 0, value);
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);

  igraph_matrix_init(&values, 0, 0);
  igraph_matrix_init(&vectors, 0, 0);
  options.mode=1;  
  igraph_sparsemat_arpack_rnsolve(&B, &options, /*storage=*/ 0,
				  &values, &vectors);
  printf("rnsolve:\n  - eigenvalues:\n    "); igraph_matrix_print(&values);
  printf("  - eigenvectors:\n    "); igraph_matrix_print(&vectors);
  igraph_matrix_destroy(&values);
  igraph_matrix_destroy(&vectors);

  igraph_vector_init(&values2, 0);
  igraph_matrix_init(&vectors, 0, 0);
  options.mode=1;  
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values2, &vectors, IGRAPH_SPARSEMAT_SOLVE_LU);
  printf("rssolve:\n  - eigenvalues:\n    "); igraph_vector_print(&values2);
  printf("  - eigenvectors:\n    "); igraph_matrix_print(&vectors);
  igraph_vector_destroy(&values2);
  igraph_matrix_destroy(&vectors);

  igraph_sparsemat_destroy(&B);
}

/* Generic test for 2x2 matrices */
void test_2x2(igraph_real_t a, igraph_real_t b, igraph_real_t c, igraph_real_t d) {
  igraph_sparsemat_t A, B;
  igraph_matrix_t values, vectors;
  igraph_vector_t values2;
  igraph_arpack_options_t options;

  igraph_arpack_options_init(&options);
  options.mode=1; options.nev=2;

  igraph_sparsemat_init(&A, 2, 2, 4);
  igraph_sparsemat_entry(&A, 0, 0, a);
  igraph_sparsemat_entry(&A, 0, 1, b);
  igraph_sparsemat_entry(&A, 1, 0, c);
  igraph_sparsemat_entry(&A, 1, 1, d);
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);

  igraph_matrix_init(&values, 0, 0);
  igraph_matrix_init(&vectors, 0, 0);
  igraph_sparsemat_arpack_rnsolve(&B, &options, /*storage=*/ 0,
				  &values, &vectors);
  printf("rnsolve:\n  - eigenvalues:\n    "); igraph_matrix_print(&values);
  printf("  - eigenvectors:\n    "); igraph_matrix_print(&vectors);
  igraph_matrix_destroy(&values);
  igraph_matrix_destroy(&vectors);

  if (b == c) {
    igraph_vector_init(&values2, 0);
    igraph_matrix_init(&vectors, 0, 0);
    igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
            &values2, &vectors, IGRAPH_SPARSEMAT_SOLVE_QR);
    printf("rssolve:\n  - eigenvalues:\n    "); igraph_vector_print(&values2);
    printf("  - eigenvectors:\n    "); igraph_matrix_print(&vectors);
    igraph_vector_destroy(&values2);
    igraph_matrix_destroy(&vectors);
  }

  igraph_sparsemat_destroy(&B);
}

int main() {

  igraph_sparsemat_t A, B;
  igraph_matrix_t vectors, values2;
  igraph_vector_t values;
  long int i;
  igraph_arpack_options_t options;
  igraph_real_t min, max;
  igraph_t g1, g2, g3;

  /***********************************************************************/

  /* Identity matrix */
#define DIM 10
  igraph_sparsemat_init(&A, DIM, DIM, DIM);
  for (i=0; i<DIM; i++) {
    igraph_sparsemat_entry(&A, i, i, 1.0);
  }
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);

  igraph_vector_init(&values, 0);
  igraph_arpack_options_init(&options);

  options.mode=1;
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0, 
				  &values, /*vectors=*/ 0, /*solvemethod=*/0);
  if (VECTOR(values)[0] != 1.0) { return 1; }

  options.mode=3;
  options.sigma=2;
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0, 
				  &values, /*vectors=*/ 0,
				  IGRAPH_SPARSEMAT_SOLVE_LU);
  if (VECTOR(values)[0] != 1.0) { return 21; }
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0, 
				  &values, /*vectors=*/ 0,
				  IGRAPH_SPARSEMAT_SOLVE_QR);
  if (VECTOR(values)[0] != 1.0) { return 31; }

  igraph_vector_destroy(&values);
  igraph_sparsemat_destroy(&B);  

#undef DIM

  /***********************************************************************/

  /* Diagonal matrix */
#define DIM 10
  igraph_sparsemat_init(&A, DIM, DIM, DIM);
  for (i=0; i<DIM; i++) {
    igraph_sparsemat_entry(&A, i, i, i+1.0);
  }
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);
  
  igraph_vector_init(&values, 0);
  igraph_matrix_init(&vectors, 0, 0);
  
  options.mode=1;
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values, /*vectors=*/ &vectors, 
				  /*solvemethod=*/ 0);
  if ( fabs(VECTOR(values)[0] - DIM) > EPS ) {
    printf("VECTOR(values)[0] numerical precision is only %g, should be %g",
			fabs((double)VECTOR(values)[0]-DIM), EPS);
	return 2;
  }

  if ( fabs(fabs(MATRIX(vectors, DIM-1, 0)) - 1.0) > EPS) { return 3; }
  MATRIX(vectors, DIM-1, 0) = 0.0;
  igraph_matrix_minmax(&vectors, &min, &max);
  if (fabs(min) > EPS) { return 3; }
  if (fabs(max) > EPS) { return 3; }

  options.mode=3;
  options.sigma=11;
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values, /*vectors=*/ &vectors,
				  IGRAPH_SPARSEMAT_SOLVE_LU);
  if ( fabs(VECTOR(values)[0] - DIM) > EPS ) {
    printf("VECTOR(values)[0] numerical precision is only %g, should be %g",
			fabs((double)VECTOR(values)[0]-DIM), EPS);
	return 22;
  }
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values, /*vectors=*/ &vectors,
				  IGRAPH_SPARSEMAT_SOLVE_QR);
  if ( fabs(VECTOR(values)[0] - DIM) > EPS ) {
    printf("VECTOR(values)[0] numerical precision is only %g, should be %g",
			fabs((double)VECTOR(values)[0]-DIM), EPS);
    return 32;
  }

  if ( fabs(fabs(MATRIX(vectors, DIM-1, 0)) - 1.0) > EPS) { return 23; }
  MATRIX(vectors, DIM-1, 0) = 0.0;
  igraph_matrix_minmax(&vectors, &min, &max);
  if (fabs(min) > EPS) { return 23; }
  if (fabs(max) > EPS) { return 23; }
  
  igraph_vector_destroy(&values);
  igraph_matrix_destroy(&vectors);
  igraph_sparsemat_destroy(&B);
#undef DIM

  /***********************************************************************/

  /* A tree, plus a ring */
#define DIM 10
  igraph_tree(&g1, DIM, /*children=*/ 2, IGRAPH_TREE_UNDIRECTED);
  igraph_ring(&g2, DIM, IGRAPH_UNDIRECTED, /*mutual=*/ 0, /*circular=*/ 1);
  igraph_union(&g3, &g1, &g2, /*edge_map1=*/ 0, /*edge_map1=*/ 0);
  igraph_destroy(&g1);
  igraph_destroy(&g2);

  igraph_get_sparsemat(&g3, &A);
  igraph_destroy(&g3);
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);
  
  igraph_vector_init(&values, 0);
  igraph_matrix_init(&vectors, 0, 0);
  
  options.mode=1;
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values, &vectors, /*solvemethod=*/ 0);

  if (MATRIX(vectors, 0, 0) < 0.0) { 
    igraph_matrix_scale(&vectors, -1.0);
  }
  
  igraph_vector_print(&values);
  igraph_matrix_print(&vectors);

  options.mode=3;
  options.sigma=VECTOR(values)[0] * 1.1;
  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values, &vectors, 
				  IGRAPH_SPARSEMAT_SOLVE_LU);

  if (MATRIX(vectors, 0, 0) < 0.0) { 
    igraph_matrix_scale(&vectors, -1.0);
  }  
  igraph_vector_print(&values);
  igraph_matrix_print(&vectors);

  igraph_sparsemat_arpack_rssolve(&B, &options, /*storage=*/ 0,
				  &values, &vectors, 
				  IGRAPH_SPARSEMAT_SOLVE_QR);
  if (MATRIX(vectors, 0, 0) < 0.0) { 
    igraph_matrix_scale(&vectors, -1.0);
  }  
  igraph_vector_print(&values);
  igraph_matrix_print(&vectors);
  
  igraph_vector_destroy(&values);
  igraph_matrix_destroy(&vectors);
  igraph_sparsemat_destroy(&B);
#undef DIM

  printf("--\n");

  /***********************************************************************/

  /* A directed tree and a directed, mutual ring */
#define DIM 10
  igraph_tree(&g1, DIM, /*children=*/ 2, IGRAPH_TREE_OUT);
  igraph_ring(&g2, DIM, IGRAPH_DIRECTED, /*mutual=*/ 1, /*circular=*/ 1);
  igraph_union(&g3, &g1, &g2, /*edge_map1=*/ 0, /*edge_map2=*/ 0);
  igraph_destroy(&g1);
  igraph_destroy(&g2);

  igraph_get_sparsemat(&g3, &A);
  igraph_destroy(&g3);
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);
  
  igraph_matrix_init(&values2, 0, 0);
  igraph_matrix_init(&vectors, 0, 0);
  
  options.mode=1;
  igraph_sparsemat_arpack_rnsolve(&B, &options, /*storage=*/ 0,
				  &values2, &vectors);

  if (MATRIX(vectors, 0, 0) < 0.0) { 
    igraph_matrix_scale(&vectors, -1.0);
  }
  
  igraph_matrix_print(&values2);
  igraph_matrix_print(&vectors);
  
  igraph_matrix_destroy(&values2);
  igraph_matrix_destroy(&vectors);
  igraph_sparsemat_destroy(&B);
#undef DIM

  /***********************************************************************/

  /* A small test graph */

  igraph_small(&g1, 11, IGRAPH_DIRECTED,
	       0, 1, 1, 3, 1, 8, 2, 10, 3, 6, 3, 10, 4, 2, 5, 4, 
	       6, 1, 6, 4, 7, 9, 8, 5, 8, 7, 9, 8, 10, 0,
	       -1);

  igraph_get_sparsemat(&g1, &A);
  igraph_destroy(&g1);
  igraph_sparsemat_compress(&A, &B);
  igraph_sparsemat_destroy(&A);
  
  igraph_matrix_init(&values2, 0, 0);
  igraph_matrix_init(&vectors, 0, 0);

  options.mode=1;  
  igraph_sparsemat_arpack_rnsolve(&B, &options, /*storage=*/ 0,
				  &values2, &vectors);

  if (MATRIX(vectors, 0, 0) < 0.0) { 
    igraph_matrix_scale(&vectors, -1.0);
  }
  
  igraph_matrix_destroy(&values2);
  igraph_matrix_destroy(&vectors);
  igraph_sparsemat_destroy(&B);

  /***********************************************************************/

  /* Testing the special case solver for 1x1 matrices */
  printf("--\n");
  test_1x1(2);
  test_1x1(0);
  test_1x1(-3);

  /***********************************************************************/

  /* Testing the special case solver for 2x2 matrices */
  printf("--\n");
  test_2x2(1, 2, 2, 4);      /* symmetric */
  test_2x2(1, 2, 3, 4);      /* non-symmetric, real eigenvalues */
  test_2x2(1, -5, 10, 4);    /* non-symmetric, complex eigenvalues */
  test_2x2(0, 0, 0, 0);      /* symmetric, pathological */

  return 0;
}