Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
amesos_colamd_l.c
Go to the documentation of this file.
1 /* ========================================================================== */
2 /* === colamd/symamd - a sparse matrix column ordering algorithm ============ */
3 /* ========================================================================== */
4 
5 /* COLAMD / SYMAMD
6 
7  colamd: an approximate minimum degree column ordering algorithm,
8  for LU factorization of symmetric or unsymmetric matrices,
9  QR factorization, least squares, interior point methods for
10  linear programming problems, and other related problems.
11 
12  symamd: an approximate minimum degree ordering algorithm for Cholesky
13  factorization of symmetric matrices.
14 
15  Purpose:
16 
17  Colamd computes a permutation Q such that the Cholesky factorization of
18  (AQ)'(AQ) has less fill-in and requires fewer floating point operations
19  than A'A. This also provides a good ordering for sparse partial
20  pivoting methods, P(AQ) = LU, where Q is computed prior to numerical
21  factorization, and P is computed during numerical factorization via
22  conventional partial pivoting with row interchanges. Colamd is the
23  column ordering method used in SuperLU, part of the ScaLAPACK library.
24  It is also available as built-in function in MATLAB Version 6,
25  available from MathWorks, Inc. (http://www.mathworks.com). This
26  routine can be used in place of colmmd in MATLAB.
27 
28  Symamd computes a permutation P of a symmetric matrix A such that the
29  Cholesky factorization of PAP' has less fill-in and requires fewer
30  floating point operations than A. Symamd constructs a matrix M such
31  that M'M has the same nonzero pattern of A, and then orders the columns
32  of M using colmmd. The column ordering of M is then returned as the
33  row and column ordering P of A.
34 
35  Authors:
36 
37  The authors of the code itself are Stefan I. Larimore and Timothy A.
38  Davis (davis at cise.ufl.edu), University of Florida. The algorithm was
39  developed in collaboration with John Gilbert, Xerox PARC, and Esmond
40  Ng, Oak Ridge National Laboratory.
41 
42  Acknowledgements:
43 
44  This work was supported by the National Science Foundation, under
45  grants DMS-9504974 and DMS-9803599.
46 
47  Copyright and License:
48 
49  Copyright (c) 1998-2007, Timothy A. Davis, All Rights Reserved.
50  COLAMD is also available under alternate licenses, contact T. Davis
51  for details.
52 
53  This library is free software; you can redistribute it and/or
54  modify it under the terms of the GNU Lesser General Public
55  License as published by the Free Software Foundation; either
56  version 2.1 of the License, or (at your option) any later version.
57 
58  This library is distributed in the hope that it will be useful,
59  but WITHOUT ANY WARRANTY; without even the implied warranty of
60  MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
61  Lesser General Public License for more details.
62 
63  You should have received a copy of the GNU Lesser General Public
64  License along with this library; if not, write to the Free Software
65  Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
66  USA
67 
68  Permission is hereby granted to use or copy this program under the
69  terms of the GNU LGPL, provided that the Copyright, this License,
70  and the Availability of the original version is retained on all copies.
71  User documentation of any code that uses this code or any modified
72  version of this code must cite the Copyright, this License, the
73  Availability note, and "Used by permission." Permission to modify
74  the code and to distribute modified code is granted, provided the
75  Copyright, this License, and the Availability note are retained,
76  and a notice that the code was modified is included.
77 
78  Availability:
79 
80  The colamd/symamd library is available at
81 
82  http://www.cise.ufl.edu/research/sparse/colamd/
83 
84  This is the http://www.cise.ufl.edu/research/sparse/colamd/colamd.c
85  file. It requires the colamd.h file. It is required by the colamdmex.c
86  and symamdmex.c files, for the MATLAB interface to colamd and symamd.
87  Appears as ACM Algorithm 836.
88 
89  See the ChangeLog file for changes since Version 1.0.
90 
91  References:
92 
93  T. A. Davis, J. R. Gilbert, S. Larimore, E. Ng, An approximate column
94  minimum degree ordering algorithm, ACM Transactions on Mathematical
95  Software, vol. 30, no. 3., pp. 353-376, 2004.
96 
97  T. A. Davis, J. R. Gilbert, S. Larimore, E. Ng, Algorithm 836: COLAMD,
98  an approximate column minimum degree ordering algorithm, ACM
99  Transactions on Mathematical Software, vol. 30, no. 3., pp. 377-380,
100  2004.
101 
102 */
103 
104 /* ========================================================================== */
105 /* === Description of user-callable routines ================================ */
106 /* ========================================================================== */
107 
108 /* COLAMD includes both int and UF_long versions of all its routines. The
109  * description below is for the int version. For UF_long, all int arguments
110  * become UF_long. UF_long is normally defined as long, except for WIN64.
111 
112  ----------------------------------------------------------------------------
113  colamd_recommended:
114  ----------------------------------------------------------------------------
115 
116  C syntax:
117 
118  #include "colamd.h"
119  size_t colamd_recommended (int nnz, int n_row, int n_col) ;
120  size_t colamd_l_recommended (UF_long nnz, UF_long n_row,
121  UF_long n_col) ;
122 
123  Purpose:
124 
125  Returns recommended value of Alen for use by colamd. Returns 0
126  if any input argument is negative. The use of this routine
127  is optional. Not needed for symamd, which dynamically allocates
128  its own memory.
129 
130  Note that in v2.4 and earlier, these routines returned int or long.
131  They now return a value of type size_t.
132 
133  Arguments (all input arguments):
134 
135  int nnz ; Number of nonzeros in the matrix A. This must
136  be the same value as p [n_col] in the call to
137  colamd - otherwise you will get a wrong value
138  of the recommended memory to use.
139 
140  int n_row ; Number of rows in the matrix A.
141 
142  int n_col ; Number of columns in the matrix A.
143 
144  ----------------------------------------------------------------------------
145  colamd_set_defaults:
146  ----------------------------------------------------------------------------
147 
148  C syntax:
149 
150  #include "colamd.h"
151  colamd_set_defaults (double knobs [COLAMD_KNOBS]) ;
152  colamd_l_set_defaults (double knobs [COLAMD_KNOBS]) ;
153 
154  Purpose:
155 
156  Sets the default parameters. The use of this routine is optional.
157 
158  Arguments:
159 
160  double knobs [COLAMD_KNOBS] ; Output only.
161 
162  NOTE: the meaning of the dense row/col knobs has changed in v2.4
163 
164  knobs [0] and knobs [1] control dense row and col detection:
165 
166  Colamd: rows with more than
167  max (16, knobs [COLAMD_DENSE_ROW] * sqrt (n_col))
168  entries are removed prior to ordering. Columns with more than
169  max (16, knobs [COLAMD_DENSE_COL] * sqrt (MIN (n_row,n_col)))
170  entries are removed prior to
171  ordering, and placed last in the output column ordering.
172 
173  Symamd: uses only knobs [COLAMD_DENSE_ROW], which is knobs [0].
174  Rows and columns with more than
175  max (16, knobs [COLAMD_DENSE_ROW] * sqrt (n))
176  entries are removed prior to ordering, and placed last in the
177  output ordering.
178 
179  COLAMD_DENSE_ROW and COLAMD_DENSE_COL are defined as 0 and 1,
180  respectively, in colamd.h. Default values of these two knobs
181  are both 10. Currently, only knobs [0] and knobs [1] are
182  used, but future versions may use more knobs. If so, they will
183  be properly set to their defaults by the future version of
184  colamd_set_defaults, so that the code that calls colamd will
185  not need to change, assuming that you either use
186  colamd_set_defaults, or pass a (double *) NULL pointer as the
187  knobs array to colamd or symamd.
188 
189  knobs [2]: aggressive absorption
190 
191  knobs [COLAMD_AGGRESSIVE] controls whether or not to do
192  aggressive absorption during the ordering. Default is TRUE.
193 
194 
195  ----------------------------------------------------------------------------
196  colamd:
197  ----------------------------------------------------------------------------
198 
199  C syntax:
200 
201  #include "colamd.h"
202  int colamd (int n_row, int n_col, int Alen, int *A, int *p,
203  double knobs [COLAMD_KNOBS], int stats [COLAMD_STATS]) ;
204  UF_long colamd_l (UF_long n_row, UF_long n_col, UF_long Alen,
205  UF_long *A, UF_long *p, double knobs [COLAMD_KNOBS],
206  UF_long stats [COLAMD_STATS]) ;
207 
208  Purpose:
209 
210  Computes a column ordering (Q) of A such that P(AQ)=LU or
211  (AQ)'AQ=LL' have less fill-in and require fewer floating point
212  operations than factorizing the unpermuted matrix A or A'A,
213  respectively.
214 
215  Returns:
216 
217  TRUE (1) if successful, FALSE (0) otherwise.
218 
219  Arguments:
220 
221  int n_row ; Input argument.
222 
223  Number of rows in the matrix A.
224  Restriction: n_row >= 0.
225  Colamd returns FALSE if n_row is negative.
226 
227  int n_col ; Input argument.
228 
229  Number of columns in the matrix A.
230  Restriction: n_col >= 0.
231  Colamd returns FALSE if n_col is negative.
232 
233  int Alen ; Input argument.
234 
235  Restriction (see note):
236  Alen >= 2*nnz + 6*(n_col+1) + 4*(n_row+1) + n_col
237  Colamd returns FALSE if these conditions are not met.
238 
239  Note: this restriction makes an modest assumption regarding
240  the size of the two typedef's structures in colamd.h.
241  We do, however, guarantee that
242 
243  Alen >= colamd_recommended (nnz, n_row, n_col)
244 
245  will be sufficient. Note: the macro version does not check
246  for integer overflow, and thus is not recommended. Use
247  the colamd_recommended routine instead.
248 
249  int A [Alen] ; Input argument, undefined on output.
250 
251  A is an integer array of size Alen. Alen must be at least as
252  large as the bare minimum value given above, but this is very
253  low, and can result in excessive run time. For best
254  performance, we recommend that Alen be greater than or equal to
255  colamd_recommended (nnz, n_row, n_col), which adds
256  nnz/5 to the bare minimum value given above.
257 
258  On input, the row indices of the entries in column c of the
259  matrix are held in A [(p [c]) ... (p [c+1]-1)]. The row indices
260  in a given column c need not be in ascending order, and
261  duplicate row indices may be be present. However, colamd will
262  work a little faster if both of these conditions are met
263  (Colamd puts the matrix into this format, if it finds that the
264  the conditions are not met).
265 
266  The matrix is 0-based. That is, rows are in the range 0 to
267  n_row-1, and columns are in the range 0 to n_col-1. Colamd
268  returns FALSE if any row index is out of range.
269 
270  The contents of A are modified during ordering, and are
271  undefined on output.
272 
273  int p [n_col+1] ; Both input and output argument.
274 
275  p is an integer array of size n_col+1. On input, it holds the
276  "pointers" for the column form of the matrix A. Column c of
277  the matrix A is held in A [(p [c]) ... (p [c+1]-1)]. The first
278  entry, p [0], must be zero, and p [c] <= p [c+1] must hold
279  for all c in the range 0 to n_col-1. The value p [n_col] is
280  thus the total number of entries in the pattern of the matrix A.
281  Colamd returns FALSE if these conditions are not met.
282 
283  On output, if colamd returns TRUE, the array p holds the column
284  permutation (Q, for P(AQ)=LU or (AQ)'(AQ)=LL'), where p [0] is
285  the first column index in the new ordering, and p [n_col-1] is
286  the last. That is, p [k] = j means that column j of A is the
287  kth pivot column, in AQ, where k is in the range 0 to n_col-1
288  (p [0] = j means that column j of A is the first column in AQ).
289 
290  If colamd returns FALSE, then no permutation is returned, and
291  p is undefined on output.
292 
293  double knobs [COLAMD_KNOBS] ; Input argument.
294 
295  See colamd_set_defaults for a description.
296 
297  int stats [COLAMD_STATS] ; Output argument.
298 
299  Statistics on the ordering, and error status.
300  See colamd.h for related definitions.
301  Colamd returns FALSE if stats is not present.
302 
303  stats [0]: number of dense or empty rows ignored.
304 
305  stats [1]: number of dense or empty columns ignored (and
306  ordered last in the output permutation p)
307  Note that a row can become "empty" if it
308  contains only "dense" and/or "empty" columns,
309  and similarly a column can become "empty" if it
310  only contains "dense" and/or "empty" rows.
311 
312  stats [2]: number of garbage collections performed.
313  This can be excessively high if Alen is close
314  to the minimum required value.
315 
316  stats [3]: status code. < 0 is an error code.
317  > 1 is a warning or notice.
318 
319  0 OK. Each column of the input matrix contained
320  row indices in increasing order, with no
321  duplicates.
322 
323  1 OK, but columns of input matrix were jumbled
324  (unsorted columns or duplicate entries). Colamd
325  had to do some extra work to sort the matrix
326  first and remove duplicate entries, but it
327  still was able to return a valid permutation
328  (return value of colamd was TRUE).
329 
330  stats [4]: highest numbered column that
331  is unsorted or has duplicate
332  entries.
333  stats [5]: last seen duplicate or
334  unsorted row index.
335  stats [6]: number of duplicate or
336  unsorted row indices.
337 
338  -1 A is a null pointer
339 
340  -2 p is a null pointer
341 
342  -3 n_row is negative
343 
344  stats [4]: n_row
345 
346  -4 n_col is negative
347 
348  stats [4]: n_col
349 
350  -5 number of nonzeros in matrix is negative
351 
352  stats [4]: number of nonzeros, p [n_col]
353 
354  -6 p [0] is nonzero
355 
356  stats [4]: p [0]
357 
358  -7 A is too small
359 
360  stats [4]: required size
361  stats [5]: actual size (Alen)
362 
363  -8 a column has a negative number of entries
364 
365  stats [4]: column with < 0 entries
366  stats [5]: number of entries in col
367 
368  -9 a row index is out of bounds
369 
370  stats [4]: column with bad row index
371  stats [5]: bad row index
372  stats [6]: n_row, # of rows of matrx
373 
374  -10 (unused; see symamd.c)
375 
376  -999 (unused; see symamd.c)
377 
378  Future versions may return more statistics in the stats array.
379 
380  Example:
381 
382  See http://www.cise.ufl.edu/research/sparse/colamd/example.c
383  for a complete example.
384 
385  To order the columns of a 5-by-4 matrix with 11 nonzero entries in
386  the following nonzero pattern
387 
388  x 0 x 0
389  x 0 x x
390  0 x x 0
391  0 0 x x
392  x x 0 0
393 
394  with default knobs and no output statistics, do the following:
395 
396  #include "colamd.h"
397  #define ALEN 100
398  int A [ALEN] = {0, 1, 4, 2, 4, 0, 1, 2, 3, 1, 3} ;
399  int p [ ] = {0, 3, 5, 9, 11} ;
400  int stats [COLAMD_STATS] ;
401  colamd (5, 4, ALEN, A, p, (double *) NULL, stats) ;
402 
403  The permutation is returned in the array p, and A is destroyed.
404 
405  ----------------------------------------------------------------------------
406  symamd:
407  ----------------------------------------------------------------------------
408 
409  C syntax:
410 
411  #include "colamd.h"
412  int symamd (int n, int *A, int *p, int *perm,
413  double knobs [COLAMD_KNOBS], int stats [COLAMD_STATS],
414  void (*allocate) (size_t, size_t), void (*release) (void *)) ;
415  UF_long symamd_l (UF_long n, UF_long *A, UF_long *p, UF_long *perm,
416  double knobs [COLAMD_KNOBS], UF_long stats [COLAMD_STATS],
417  void (*allocate) (size_t, size_t), void (*release) (void *)) ;
418 
419  Purpose:
420 
421  The symamd routine computes an ordering P of a symmetric sparse
422  matrix A such that the Cholesky factorization PAP' = LL' remains
423  sparse. It is based on a column ordering of a matrix M constructed
424  so that the nonzero pattern of M'M is the same as A. The matrix A
425  is assumed to be symmetric; only the strictly lower triangular part
426  is accessed. You must pass your selected memory allocator (usually
427  calloc/free or mxCalloc/mxFree) to symamd, for it to allocate
428  memory for the temporary matrix M.
429 
430  Returns:
431 
432  TRUE (1) if successful, FALSE (0) otherwise.
433 
434  Arguments:
435 
436  int n ; Input argument.
437 
438  Number of rows and columns in the symmetrix matrix A.
439  Restriction: n >= 0.
440  Symamd returns FALSE if n is negative.
441 
442  int A [nnz] ; Input argument.
443 
444  A is an integer array of size nnz, where nnz = p [n].
445 
446  The row indices of the entries in column c of the matrix are
447  held in A [(p [c]) ... (p [c+1]-1)]. The row indices in a
448  given column c need not be in ascending order, and duplicate
449  row indices may be present. However, symamd will run faster
450  if the columns are in sorted order with no duplicate entries.
451 
452  The matrix is 0-based. That is, rows are in the range 0 to
453  n-1, and columns are in the range 0 to n-1. Symamd
454  returns FALSE if any row index is out of range.
455 
456  The contents of A are not modified.
457 
458  int p [n+1] ; Input argument.
459 
460  p is an integer array of size n+1. On input, it holds the
461  "pointers" for the column form of the matrix A. Column c of
462  the matrix A is held in A [(p [c]) ... (p [c+1]-1)]. The first
463  entry, p [0], must be zero, and p [c] <= p [c+1] must hold
464  for all c in the range 0 to n-1. The value p [n] is
465  thus the total number of entries in the pattern of the matrix A.
466  Symamd returns FALSE if these conditions are not met.
467 
468  The contents of p are not modified.
469 
470  int perm [n+1] ; Output argument.
471 
472  On output, if symamd returns TRUE, the array perm holds the
473  permutation P, where perm [0] is the first index in the new
474  ordering, and perm [n-1] is the last. That is, perm [k] = j
475  means that row and column j of A is the kth column in PAP',
476  where k is in the range 0 to n-1 (perm [0] = j means
477  that row and column j of A are the first row and column in
478  PAP'). The array is used as a workspace during the ordering,
479  which is why it must be of length n+1, not just n.
480 
481  double knobs [COLAMD_KNOBS] ; Input argument.
482 
483  See colamd_set_defaults for a description.
484 
485  int stats [COLAMD_STATS] ; Output argument.
486 
487  Statistics on the ordering, and error status.
488  See colamd.h for related definitions.
489  Symamd returns FALSE if stats is not present.
490 
491  stats [0]: number of dense or empty row and columns ignored
492  (and ordered last in the output permutation
493  perm). Note that a row/column can become
494  "empty" if it contains only "dense" and/or
495  "empty" columns/rows.
496 
497  stats [1]: (same as stats [0])
498 
499  stats [2]: number of garbage collections performed.
500 
501  stats [3]: status code. < 0 is an error code.
502  > 1 is a warning or notice.
503 
504  0 OK. Each column of the input matrix contained
505  row indices in increasing order, with no
506  duplicates.
507 
508  1 OK, but columns of input matrix were jumbled
509  (unsorted columns or duplicate entries). Symamd
510  had to do some extra work to sort the matrix
511  first and remove duplicate entries, but it
512  still was able to return a valid permutation
513  (return value of symamd was TRUE).
514 
515  stats [4]: highest numbered column that
516  is unsorted or has duplicate
517  entries.
518  stats [5]: last seen duplicate or
519  unsorted row index.
520  stats [6]: number of duplicate or
521  unsorted row indices.
522 
523  -1 A is a null pointer
524 
525  -2 p is a null pointer
526 
527  -3 (unused, see colamd.c)
528 
529  -4 n is negative
530 
531  stats [4]: n
532 
533  -5 number of nonzeros in matrix is negative
534 
535  stats [4]: # of nonzeros (p [n]).
536 
537  -6 p [0] is nonzero
538 
539  stats [4]: p [0]
540 
541  -7 (unused)
542 
543  -8 a column has a negative number of entries
544 
545  stats [4]: column with < 0 entries
546  stats [5]: number of entries in col
547 
548  -9 a row index is out of bounds
549 
550  stats [4]: column with bad row index
551  stats [5]: bad row index
552  stats [6]: n_row, # of rows of matrx
553 
554  -10 out of memory (unable to allocate temporary
555  workspace for M or count arrays using the
556  "allocate" routine passed into symamd).
557 
558  Future versions may return more statistics in the stats array.
559 
560  void * (*allocate) (size_t, size_t)
561 
562  A pointer to a function providing memory allocation. The
563  allocated memory must be returned initialized to zero. For a
564  C application, this argument should normally be a pointer to
565  calloc. For a MATLAB mexFunction, the routine mxCalloc is
566  passed instead.
567 
568  void (*release) (size_t, size_t)
569 
570  A pointer to a function that frees memory allocated by the
571  memory allocation routine above. For a C application, this
572  argument should normally be a pointer to free. For a MATLAB
573  mexFunction, the routine mxFree is passed instead.
574 
575 
576  ----------------------------------------------------------------------------
577  colamd_report:
578  ----------------------------------------------------------------------------
579 
580  C syntax:
581 
582  #include "colamd.h"
583  colamd_report (int stats [COLAMD_STATS]) ;
584  colamd_l_report (UF_long stats [COLAMD_STATS]) ;
585 
586  Purpose:
587 
588  Prints the error status and statistics recorded in the stats
589  array on the standard error output (for a standard C routine)
590  or on the MATLAB output (for a mexFunction).
591 
592  Arguments:
593 
594  int stats [COLAMD_STATS] ; Input only. Statistics from colamd.
595 
596 
597  ----------------------------------------------------------------------------
598  symamd_report:
599  ----------------------------------------------------------------------------
600 
601  C syntax:
602 
603  #include "colamd.h"
604  symamd_report (int stats [COLAMD_STATS]) ;
605  symamd_l_report (UF_long stats [COLAMD_STATS]) ;
606 
607  Purpose:
608 
609  Prints the error status and statistics recorded in the stats
610  array on the standard error output (for a standard C routine)
611  or on the MATLAB output (for a mexFunction).
612 
613  Arguments:
614 
615  int stats [COLAMD_STATS] ; Input only. Statistics from symamd.
616 
617 
618 */
619 
620 /* ========================================================================== */
621 /* === Scaffolding code definitions ======================================== */
622 /* ========================================================================== */
623 
624 /* Ensure that debugging is turned off: */
625 #ifndef NDEBUG
626 #define NDEBUG
627 #endif
628 
629 /* turn on debugging by uncommenting the following line
630  #undef NDEBUG
631 */
632 
633 /*
634  Our "scaffolding code" philosophy: In our opinion, well-written library
635  code should keep its "debugging" code, and just normally have it turned off
636  by the compiler so as not to interfere with performance. This serves
637  several purposes:
638 
639  (1) assertions act as comments to the reader, telling you what the code
640  expects at that point. All assertions will always be true (unless
641  there really is a bug, of course).
642 
643  (2) leaving in the scaffolding code assists anyone who would like to modify
644  the code, or understand the algorithm (by reading the debugging output,
645  one can get a glimpse into what the code is doing).
646 
647  (3) (gasp!) for actually finding bugs. This code has been heavily tested
648  and "should" be fully functional and bug-free ... but you never know...
649 
650  The code will become outrageously slow when debugging is
651  enabled. To control the level of debugging output, set an environment
652  variable D to 0 (little), 1 (some), 2, 3, or 4 (lots). When debugging,
653  you should see the following message on the standard output:
654 
655  colamd: debug version, D = 1 (THIS WILL BE SLOW!)
656 
657  or a similar message for symamd. If you don't, then debugging has not
658  been enabled.
659 
660 */
661 
662 /* ========================================================================== */
663 /* === Include files ======================================================== */
664 /* ========================================================================== */
665 
666 /* This file should make the long int version of COLAMD */
667 #define DLONG 1
668 
669 #include "amesos_colamd.h"
670 #include <limits.h>
671 #include <math.h>
672 
673 #ifdef MATLAB_MEX_FILE
674 #include "mex.h"
675 #include "matrix.h"
676 #endif /* MATLAB_MEX_FILE */
677 
678 #if !defined (NPRINT) || !defined (NDEBUG)
679 #include <stdio.h>
680 #endif
681 
682 #ifndef NULL
683 #define NULL ((void *) 0)
684 #endif
685 
686 /* ========================================================================== */
687 /* === int or UF_long ======================================================= */
688 /* ========================================================================== */
689 
690 /* define UF_long */
691 #include "amesos_UFconfig.h"
692 
693 #ifdef DLONG
694 
695 #define Int UF_long
696 #define ID UF_long_id
697 #define Int_MAX UF_long_max
698 
699 #define COLAMD_recommended amesos_colamd_l_recommended
700 #define COLAMD_set_defaults amesos_colamd_l_set_defaults
701 #define COLAMD_MAIN amesos_colamd_l
702 #define SYMAMD_MAIN amesos_symamd_l
703 #define COLAMD_report amesos_colamd_l_report
704 #define SYMAMD_report amesos_symamd_l_report
705 
706 #else
707 
708 #define Int int
709 #define ID "%d"
710 #define Int_MAX INT_MAX
711 
712 #define COLAMD_recommended amesos_colamd_recommended
713 #define COLAMD_set_defaults amesos_colamd_set_defaults
714 #define COLAMD_MAIN amesos_colamd
715 #define SYMAMD_MAIN amesos_symamd
716 #define COLAMD_report amesos_colamd_report
717 #define SYMAMD_report amesos_symamd_report
718 
719 #endif
720 
721 /* ========================================================================== */
722 /* === Row and Column structures ============================================ */
723 /* ========================================================================== */
724 
725 /* User code that makes use of the colamd/symamd routines need not directly */
726 /* reference these structures. They are used only for colamd_recommended. */
727 
728 typedef struct Colamd_Col_struct
729 {
730  Int start ; /* index for A of first row in this column, or DEAD */
731  /* if column is dead */
732  Int length ; /* number of rows in this column */
733  union
734  {
735  Int thickness ; /* number of original columns represented by this */
736  /* col, if the column is alive */
737  Int parent ; /* parent in parent tree super-column structure, if */
738  /* the column is dead */
739  } shared1 ;
740  union
741  {
742  Int score ; /* the score used to maintain heap, if col is alive */
743  Int order ; /* pivot ordering of this column, if col is dead */
744  } shared2 ;
745  union
746  {
747  Int headhash ; /* head of a hash bucket, if col is at the head of */
748  /* a degree list */
749  Int hash ; /* hash value, if col is not in a degree list */
750  Int prev ; /* previous column in degree list, if col is in a */
751  /* degree list (but not at the head of a degree list) */
752  } shared3 ;
753  union
754  {
755  Int degree_next ; /* next column, if col is in a degree list */
756  Int hash_next ; /* next column, if col is in a hash list */
757  } shared4 ;
758 
759 } Colamd_Col ;
760 
761 typedef struct Colamd_Row_struct
762 {
763  Int start ; /* index for A of first col in this row */
764  Int length ; /* number of principal columns in this row */
765  union
766  {
767  Int degree ; /* number of principal & non-principal columns in row */
768  Int p ; /* used as a row pointer in init_rows_cols () */
769  } shared1 ;
770  union
771  {
772  Int mark ; /* for computing set differences and marking dead rows*/
773  Int first_column ;/* first column in row (used in garbage collection) */
774  } shared2 ;
775 
776 } Colamd_Row ;
777 
778 /* ========================================================================== */
779 /* === Definitions ========================================================== */
780 /* ========================================================================== */
781 
782 /* Routines are either PUBLIC (user-callable) or PRIVATE (not user-callable) */
783 #define PUBLIC
784 #define PRIVATE static
785 
786 #define DENSE_DEGREE(alpha,n) \
787  ((Int) MAX (16.0, (alpha) * sqrt ((double) (n))))
788 
789 #define MAX(a,b) (((a) > (b)) ? (a) : (b))
790 #define MIN(a,b) (((a) < (b)) ? (a) : (b))
791 
792 #define ONES_COMPLEMENT(r) (-(r)-1)
793 
794 /* -------------------------------------------------------------------------- */
795 /* Change for version 2.1: define TRUE and FALSE only if not yet defined */
796 /* -------------------------------------------------------------------------- */
797 
798 #ifndef TRUE
799 #define TRUE (1)
800 #endif
801 
802 #ifndef FALSE
803 #define FALSE (0)
804 #endif
805 
806 /* -------------------------------------------------------------------------- */
807 
808 #define EMPTY (-1)
809 
810 /* Row and column status */
811 #define ALIVE (0)
812 #define DEAD (-1)
813 
814 /* Column status */
815 #define DEAD_PRINCIPAL (-1)
816 #define DEAD_NON_PRINCIPAL (-2)
817 
818 /* Macros for row and column status update and checking. */
819 #define ROW_IS_DEAD(r) ROW_IS_MARKED_DEAD (Row[r].shared2.mark)
820 #define ROW_IS_MARKED_DEAD(row_mark) (row_mark < ALIVE)
821 #define ROW_IS_ALIVE(r) (Row [r].shared2.mark >= ALIVE)
822 #define COL_IS_DEAD(c) (Col [c].start < ALIVE)
823 #define COL_IS_ALIVE(c) (Col [c].start >= ALIVE)
824 #define COL_IS_DEAD_PRINCIPAL(c) (Col [c].start == DEAD_PRINCIPAL)
825 #define KILL_ROW(r) { Row [r].shared2.mark = DEAD ; }
826 #define KILL_PRINCIPAL_COL(c) { Col [c].start = DEAD_PRINCIPAL ; }
827 #define KILL_NON_PRINCIPAL_COL(c) { Col [c].start = DEAD_NON_PRINCIPAL ; }
828 
829 /* ========================================================================== */
830 /* === Colamd reporting mechanism =========================================== */
831 /* ========================================================================== */
832 
833 #if defined (MATLAB_MEX_FILE) || defined (MATHWORKS)
834 /* In MATLAB, matrices are 1-based to the user, but 0-based internally */
835 #define INDEX(i) ((i)+1)
836 #else
837 /* In C, matrices are 0-based and indices are reported as such in *_report */
838 #define INDEX(i) (i)
839 #endif
840 
841 /* All output goes through the PRINTF macro. */
842 #define PRINTF(params) { if (amesos_colamd_printf != NULL) (void) amesos_colamd_printf params ; }
843 
844 /* ========================================================================== */
845 /* === Prototypes of PRIVATE routines ======================================= */
846 /* ========================================================================== */
847 
849 (
850  Int n_row,
851  Int n_col,
852  Colamd_Row Row [],
853  Colamd_Col Col [],
854  Int A [],
855  Int p [],
856  Int stats [COLAMD_STATS]
857 ) ;
858 
860 (
861  Int n_row,
862  Int n_col,
863  Colamd_Row Row [],
864  Colamd_Col Col [],
865  Int A [],
866  Int head [],
867  double knobs [COLAMD_KNOBS],
868  Int *p_n_row2,
869  Int *p_n_col2,
870  Int *p_max_deg
871 ) ;
872 
874 (
875  Int n_row,
876  Int n_col,
877  Int Alen,
878  Colamd_Row Row [],
879  Colamd_Col Col [],
880  Int A [],
881  Int head [],
882  Int n_col2,
883  Int max_deg,
884  Int pfree,
885  Int aggressive
886 ) ;
887 
889 (
890  Int n_col,
891  Colamd_Col Col [],
892  Int p []
893 ) ;
894 
896 (
897 
898 #ifndef NDEBUG
899  Int n_col,
900  Colamd_Row Row [],
901 #endif /* NDEBUG */
902 
903  Colamd_Col Col [],
904  Int A [],
905  Int head [],
906  Int row_start,
907  Int row_length
908 ) ;
909 
911 (
912  Int n_row,
913  Int n_col,
914  Colamd_Row Row [],
915  Colamd_Col Col [],
916  Int A [],
917  Int *pfree
918 ) ;
919 
921 (
922  Int tag_mark,
923  Int max_mark,
924  Int n_row,
925  Colamd_Row Row []
926 ) ;
927 
929 (
930  char *method,
931  Int stats [COLAMD_STATS]
932 ) ;
933 
934 /* ========================================================================== */
935 /* === Debugging prototypes and definitions ================================= */
936 /* ========================================================================== */
937 
938 #ifndef NDEBUG
939 
940 #include <assert.h>
941 
942 /* colamd_debug is the *ONLY* global variable, and is only */
943 /* present when debugging */
944 
945 PRIVATE Int colamd_debug = 0 ; /* debug print level */
946 
947 #define DEBUG0(params) { PRINTF (params) ; }
948 #define DEBUG1(params) { if (colamd_debug >= 1) PRINTF (params) ; }
949 #define DEBUG2(params) { if (colamd_debug >= 2) PRINTF (params) ; }
950 #define DEBUG3(params) { if (colamd_debug >= 3) PRINTF (params) ; }
951 #define DEBUG4(params) { if (colamd_debug >= 4) PRINTF (params) ; }
952 
953 #ifdef MATLAB_MEX_FILE
954 #define ASSERT(expression) (mxAssert ((expression), ""))
955 #else
956 #define ASSERT(expression) (assert (expression))
957 #endif /* MATLAB_MEX_FILE */
958 
959 PRIVATE void colamd_get_debug /* gets the debug print level from getenv */
960 (
961  char *method
962 ) ;
963 
964 PRIVATE void debug_deg_lists
965 (
966  Int n_row,
967  Int n_col,
968  Colamd_Row Row [],
969  Colamd_Col Col [],
970  Int head [],
971  Int min_score,
972  Int should,
973  Int max_deg
974 ) ;
975 
976 PRIVATE void debug_mark
977 (
978  Int n_row,
979  Colamd_Row Row [],
980  Int tag_mark,
981  Int max_mark
982 ) ;
983 
984 PRIVATE void debug_matrix
985 (
986  Int n_row,
987  Int n_col,
988  Colamd_Row Row [],
989  Colamd_Col Col [],
990  Int A []
991 ) ;
992 
993 PRIVATE void debug_structures
994 (
995  Int n_row,
996  Int n_col,
997  Colamd_Row Row [],
998  Colamd_Col Col [],
999  Int A [],
1000  Int n_col2
1001 ) ;
1002 
1003 #else /* NDEBUG */
1004 
1005 /* === No debugging ========================================================= */
1006 
1007 #define DEBUG0(params) ;
1008 #define DEBUG1(params) ;
1009 #define DEBUG2(params) ;
1010 #define DEBUG3(params) ;
1011 #define DEBUG4(params) ;
1012 
1013 #define ASSERT(expression)
1014 
1015 #endif /* NDEBUG */
1016 
1017 /* ========================================================================== */
1018 /* === USER-CALLABLE ROUTINES: ============================================== */
1019 /* ========================================================================== */
1020 
1021 /* ========================================================================== */
1022 /* === colamd_recommended =================================================== */
1023 /* ========================================================================== */
1024 
1025 /*
1026  The colamd_recommended routine returns the suggested size for Alen. This
1027  value has been determined to provide good balance between the number of
1028  garbage collections and the memory requirements for colamd. If any
1029  argument is negative, or if integer overflow occurs, a 0 is returned as an
1030  error condition. 2*nnz space is required for the row and column
1031  indices of the matrix. COLAMD_C (n_col) + COLAMD_R (n_row) space is
1032  required for the Col and Row arrays, respectively, which are internal to
1033  colamd (roughly 6*n_col + 4*n_row). An additional n_col space is the
1034  minimal amount of "elbow room", and nnz/5 more space is recommended for
1035  run time efficiency.
1036 
1037  Alen is approximately 2.2*nnz + 7*n_col + 4*n_row + 10.
1038 
1039  This function is not needed when using symamd.
1040 */
1041 
1042 /* add two values of type size_t, and check for integer overflow */
1043 static size_t t_add (size_t a, size_t b, int *ok)
1044 {
1045  (*ok) = (*ok) && ((a + b) >= MAX (a,b)) ;
1046  return ((*ok) ? (a + b) : 0) ;
1047 }
1048 
1049 /* compute a*k where k is a small integer, and check for integer overflow */
1050 static size_t t_mult (size_t a, size_t k, int *ok)
1051 {
1052  size_t i, s = 0 ;
1053  for (i = 0 ; i < k ; i++)
1054  {
1055  s = t_add (s, a, ok) ;
1056  }
1057  return (s) ;
1058 }
1059 
1060 /* size of the Col and Row structures */
1061 #define COLAMD_C(n_col,ok) \
1062  ((t_mult (t_add (n_col, 1, ok), sizeof (Colamd_Col), ok) / sizeof (Int)))
1063 
1064 #define COLAMD_R(n_row,ok) \
1065  ((t_mult (t_add (n_row, 1, ok), sizeof (Colamd_Row), ok) / sizeof (Int)))
1066 
1067 
1068 PUBLIC size_t COLAMD_recommended /* returns recommended value of Alen. */
1070  /* === Parameters ======================================================= */
1071 
1072  Int nnz, /* number of nonzeros in A */
1073  Int n_row, /* number of rows in A */
1074  Int n_col /* number of columns in A */
1075 )
1076 {
1077  size_t s, c, r ;
1078  int ok = TRUE ;
1079  if (nnz < 0 || n_row < 0 || n_col < 0)
1080  {
1081  return (0) ;
1082  }
1083  s = t_mult (nnz, 2, &ok) ; /* 2*nnz */
1084  c = COLAMD_C (n_col, &ok) ; /* size of column structures */
1085  r = COLAMD_R (n_row, &ok) ; /* size of row structures */
1086  s = t_add (s, c, &ok) ;
1087  s = t_add (s, r, &ok) ;
1088  s = t_add (s, n_col, &ok) ; /* elbow room */
1089  s = t_add (s, nnz/5, &ok) ; /* elbow room */
1090  ok = ok && (s < Int_MAX) ;
1091  return (ok ? s : 0) ;
1092 }
1093 
1094 
1095 /* ========================================================================== */
1096 /* === colamd_set_defaults ================================================== */
1097 /* ========================================================================== */
1098 
1099 /*
1100  The colamd_set_defaults routine sets the default values of the user-
1101  controllable parameters for colamd and symamd:
1102 
1103  Colamd: rows with more than max (16, knobs [0] * sqrt (n_col))
1104  entries are removed prior to ordering. Columns with more than
1105  max (16, knobs [1] * sqrt (MIN (n_row,n_col))) entries are removed
1106  prior to ordering, and placed last in the output column ordering.
1107 
1108  Symamd: Rows and columns with more than max (16, knobs [0] * sqrt (n))
1109  entries are removed prior to ordering, and placed last in the
1110  output ordering.
1111 
1112  knobs [0] dense row control
1113 
1114  knobs [1] dense column control
1115 
1116  knobs [2] if nonzero, do aggresive absorption
1117 
1118  knobs [3..19] unused, but future versions might use this
1119 
1120 */
1121 
1124  /* === Parameters ======================================================= */
1125 
1126  double knobs [COLAMD_KNOBS] /* knob array */
1127 )
1128 {
1129  /* === Local variables ================================================== */
1130 
1131  Int i ;
1132 
1133  if (!knobs)
1134  {
1135  return ; /* no knobs to initialize */
1136  }
1137  for (i = 0 ; i < COLAMD_KNOBS ; i++)
1138  {
1139  knobs [i] = 0 ;
1140  }
1141  knobs [COLAMD_DENSE_ROW] = 10 ;
1142  knobs [COLAMD_DENSE_COL] = 10 ;
1143  knobs [COLAMD_AGGRESSIVE] = TRUE ; /* default: do aggressive absorption*/
1144 }
1145 
1146 
1147 /* ========================================================================== */
1148 /* === symamd =============================================================== */
1149 /* ========================================================================== */
1150 
1151 PUBLIC Int SYMAMD_MAIN /* return TRUE if OK, FALSE otherwise */
1153  /* === Parameters ======================================================= */
1154 
1155  Int n, /* number of rows and columns of A */
1156  Int A [], /* row indices of A */
1157  Int p [], /* column pointers of A */
1158  Int perm [], /* output permutation, size n+1 */
1159  double knobs [COLAMD_KNOBS], /* parameters (uses defaults if NULL) */
1160  Int stats [COLAMD_STATS], /* output statistics and error codes */
1161  void * (*allocate) (size_t, size_t),
1162  /* pointer to calloc (ANSI C) or */
1163  /* mxCalloc (for MATLAB mexFunction) */
1164  void (*release) (void *)
1165  /* pointer to free (ANSI C) or */
1166  /* mxFree (for MATLAB mexFunction) */
1167 )
1168 {
1169  /* === Local variables ================================================== */
1170 
1171  Int *count ; /* length of each column of M, and col pointer*/
1172  Int *mark ; /* mark array for finding duplicate entries */
1173  Int *M ; /* row indices of matrix M */
1174  size_t Mlen ; /* length of M */
1175  Int n_row ; /* number of rows in M */
1176  Int nnz ; /* number of entries in A */
1177  Int i ; /* row index of A */
1178  Int j ; /* column index of A */
1179  Int k ; /* row index of M */
1180  Int mnz ; /* number of nonzeros in M */
1181  Int pp ; /* index into a column of A */
1182  Int last_row ; /* last row seen in the current column */
1183  Int length ; /* number of nonzeros in a column */
1184 
1185  double cknobs [COLAMD_KNOBS] ; /* knobs for colamd */
1186  double default_knobs [COLAMD_KNOBS] ; /* default knobs for colamd */
1187 
1188 #ifndef NDEBUG
1189  colamd_get_debug ("symamd") ;
1190 #endif /* NDEBUG */
1191 
1192  /* === Check the input arguments ======================================== */
1193 
1194  if (!stats)
1195  {
1196  DEBUG0 (("symamd: stats not present\n")) ;
1197  return (FALSE) ;
1198  }
1199  for (i = 0 ; i < COLAMD_STATS ; i++)
1200  {
1201  stats [i] = 0 ;
1202  }
1203  stats [COLAMD_STATUS] = COLAMD_OK ;
1204  stats [COLAMD_INFO1] = -1 ;
1205  stats [COLAMD_INFO2] = -1 ;
1206 
1207  if (!A)
1208  {
1210  DEBUG0 (("symamd: A not present\n")) ;
1211  return (FALSE) ;
1212  }
1213 
1214  if (!p) /* p is not present */
1215  {
1217  DEBUG0 (("symamd: p not present\n")) ;
1218  return (FALSE) ;
1219  }
1220 
1221  if (n < 0) /* n must be >= 0 */
1222  {
1224  stats [COLAMD_INFO1] = n ;
1225  DEBUG0 (("symamd: n negative %d\n", n)) ;
1226  return (FALSE) ;
1227  }
1228 
1229  nnz = p [n] ;
1230  if (nnz < 0) /* nnz must be >= 0 */
1231  {
1233  stats [COLAMD_INFO1] = nnz ;
1234  DEBUG0 (("symamd: number of entries negative %d\n", nnz)) ;
1235  return (FALSE) ;
1236  }
1237 
1238  if (p [0] != 0)
1239  {
1241  stats [COLAMD_INFO1] = p [0] ;
1242  DEBUG0 (("symamd: p[0] not zero %d\n", p [0])) ;
1243  return (FALSE) ;
1244  }
1245 
1246  /* === If no knobs, set default knobs =================================== */
1247 
1248  if (!knobs)
1249  {
1250  COLAMD_set_defaults (default_knobs) ;
1251  knobs = default_knobs ;
1252  }
1253 
1254  /* === Allocate count and mark ========================================== */
1255 
1256  count = (Int *) ((*allocate) (n+1, sizeof (Int))) ;
1257  if (!count)
1258  {
1260  DEBUG0 (("symamd: allocate count (size %d) failed\n", n+1)) ;
1261  return (FALSE) ;
1262  }
1263 
1264  mark = (Int *) ((*allocate) (n+1, sizeof (Int))) ;
1265  if (!mark)
1266  {
1268  (*release) ((void *) count) ;
1269  DEBUG0 (("symamd: allocate mark (size %d) failed\n", n+1)) ;
1270  return (FALSE) ;
1271  }
1272 
1273  /* === Compute column counts of M, check if A is valid ================== */
1274 
1275  stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
1276 
1277  for (i = 0 ; i < n ; i++)
1278  {
1279  mark [i] = -1 ;
1280  }
1281 
1282  for (j = 0 ; j < n ; j++)
1283  {
1284  last_row = -1 ;
1285 
1286  length = p [j+1] - p [j] ;
1287  if (length < 0)
1288  {
1289  /* column pointers must be non-decreasing */
1291  stats [COLAMD_INFO1] = j ;
1292  stats [COLAMD_INFO2] = length ;
1293  (*release) ((void *) count) ;
1294  (*release) ((void *) mark) ;
1295  DEBUG0 (("symamd: col %d negative length %d\n", j, length)) ;
1296  return (FALSE) ;
1297  }
1298 
1299  for (pp = p [j] ; pp < p [j+1] ; pp++)
1300  {
1301  i = A [pp] ;
1302  if (i < 0 || i >= n)
1303  {
1304  /* row index i, in column j, is out of bounds */
1306  stats [COLAMD_INFO1] = j ;
1307  stats [COLAMD_INFO2] = i ;
1308  stats [COLAMD_INFO3] = n ;
1309  (*release) ((void *) count) ;
1310  (*release) ((void *) mark) ;
1311  DEBUG0 (("symamd: row %d col %d out of bounds\n", i, j)) ;
1312  return (FALSE) ;
1313  }
1314 
1315  if (i <= last_row || mark [i] == j)
1316  {
1317  /* row index is unsorted or repeated (or both), thus col */
1318  /* is jumbled. This is a notice, not an error condition. */
1320  stats [COLAMD_INFO1] = j ;
1321  stats [COLAMD_INFO2] = i ;
1322  (stats [COLAMD_INFO3]) ++ ;
1323  DEBUG1 (("symamd: row %d col %d unsorted/duplicate\n", i, j)) ;
1324  }
1325 
1326  if (i > j && mark [i] != j)
1327  {
1328  /* row k of M will contain column indices i and j */
1329  count [i]++ ;
1330  count [j]++ ;
1331  }
1332 
1333  /* mark the row as having been seen in this column */
1334  mark [i] = j ;
1335 
1336  last_row = i ;
1337  }
1338  }
1339 
1340  /* v2.4: removed free(mark) */
1341 
1342  /* === Compute column pointers of M ===================================== */
1343 
1344  /* use output permutation, perm, for column pointers of M */
1345  perm [0] = 0 ;
1346  for (j = 1 ; j <= n ; j++)
1347  {
1348  perm [j] = perm [j-1] + count [j-1] ;
1349  }
1350  for (j = 0 ; j < n ; j++)
1351  {
1352  count [j] = perm [j] ;
1353  }
1354 
1355  /* === Construct M ====================================================== */
1356 
1357  mnz = perm [n] ;
1358  n_row = mnz / 2 ;
1359  Mlen = COLAMD_recommended (mnz, n_row, n) ;
1360  M = (Int *) ((*allocate) (Mlen, sizeof (Int))) ;
1361  DEBUG0 (("symamd: M is %d-by-%d with %d entries, Mlen = %g\n",
1362  n_row, n, mnz, (double) Mlen)) ;
1363 
1364  if (!M)
1365  {
1367  (*release) ((void *) count) ;
1368  (*release) ((void *) mark) ;
1369  DEBUG0 (("symamd: allocate M (size %g) failed\n", (double) Mlen)) ;
1370  return (FALSE) ;
1371  }
1372 
1373  k = 0 ;
1374 
1375  if (stats [COLAMD_STATUS] == COLAMD_OK)
1376  {
1377  /* Matrix is OK */
1378  for (j = 0 ; j < n ; j++)
1379  {
1380  ASSERT (p [j+1] - p [j] >= 0) ;
1381  for (pp = p [j] ; pp < p [j+1] ; pp++)
1382  {
1383  i = A [pp] ;
1384  ASSERT (i >= 0 && i < n) ;
1385  if (i > j)
1386  {
1387  /* row k of M contains column indices i and j */
1388  M [count [i]++] = k ;
1389  M [count [j]++] = k ;
1390  k++ ;
1391  }
1392  }
1393  }
1394  }
1395  else
1396  {
1397  /* Matrix is jumbled. Do not add duplicates to M. Unsorted cols OK. */
1398  DEBUG0 (("symamd: Duplicates in A.\n")) ;
1399  for (i = 0 ; i < n ; i++)
1400  {
1401  mark [i] = -1 ;
1402  }
1403  for (j = 0 ; j < n ; j++)
1404  {
1405  ASSERT (p [j+1] - p [j] >= 0) ;
1406  for (pp = p [j] ; pp < p [j+1] ; pp++)
1407  {
1408  i = A [pp] ;
1409  ASSERT (i >= 0 && i < n) ;
1410  if (i > j && mark [i] != j)
1411  {
1412  /* row k of M contains column indices i and j */
1413  M [count [i]++] = k ;
1414  M [count [j]++] = k ;
1415  k++ ;
1416  mark [i] = j ;
1417  }
1418  }
1419  }
1420  /* v2.4: free(mark) moved below */
1421  }
1422 
1423  /* count and mark no longer needed */
1424  (*release) ((void *) count) ;
1425  (*release) ((void *) mark) ; /* v2.4: free (mark) moved here */
1426  ASSERT (k == n_row) ;
1427 
1428  /* === Adjust the knobs for M =========================================== */
1429 
1430  for (i = 0 ; i < COLAMD_KNOBS ; i++)
1431  {
1432  cknobs [i] = knobs [i] ;
1433  }
1434 
1435  /* there are no dense rows in M */
1436  cknobs [COLAMD_DENSE_ROW] = -1 ;
1437  cknobs [COLAMD_DENSE_COL] = knobs [COLAMD_DENSE_ROW] ;
1438 
1439  /* === Order the columns of M =========================================== */
1440 
1441  /* v2.4: colamd cannot fail here, so the error check is removed */
1442  (void) COLAMD_MAIN (n_row, n, (Int) Mlen, M, perm, cknobs, stats) ;
1443 
1444  /* Note that the output permutation is now in perm */
1445 
1446  /* === get the statistics for symamd from colamd ======================== */
1447 
1448  /* a dense column in colamd means a dense row and col in symamd */
1449  stats [COLAMD_DENSE_ROW] = stats [COLAMD_DENSE_COL] ;
1450 
1451  /* === Free M =========================================================== */
1452 
1453  (*release) ((void *) M) ;
1454  DEBUG0 (("symamd: done.\n")) ;
1455  return (TRUE) ;
1456 
1457 }
1458 
1459 /* ========================================================================== */
1460 /* === colamd =============================================================== */
1461 /* ========================================================================== */
1462 
1463 /*
1464  The colamd routine computes a column ordering Q of a sparse matrix
1465  A such that the LU factorization P(AQ) = LU remains sparse, where P is
1466  selected via partial pivoting. The routine can also be viewed as
1467  providing a permutation Q such that the Cholesky factorization
1468  (AQ)'(AQ) = LL' remains sparse.
1469 */
1470 
1471 PUBLIC Int COLAMD_MAIN /* returns TRUE if successful, FALSE otherwise*/
1473  /* === Parameters ======================================================= */
1474 
1475  Int n_row, /* number of rows in A */
1476  Int n_col, /* number of columns in A */
1477  Int Alen, /* length of A */
1478  Int A [], /* row indices of A */
1479  Int p [], /* pointers to columns in A */
1480  double knobs [COLAMD_KNOBS],/* parameters (uses defaults if NULL) */
1481  Int stats [COLAMD_STATS] /* output statistics and error codes */
1482 )
1483 {
1484  /* === Local variables ================================================== */
1485 
1486  Int i ; /* loop index */
1487  Int nnz ; /* nonzeros in A */
1488  size_t Row_size ; /* size of Row [], in integers */
1489  size_t Col_size ; /* size of Col [], in integers */
1490  size_t need ; /* minimum required length of A */
1491  Colamd_Row *Row ; /* pointer into A of Row [0..n_row] array */
1492  Colamd_Col *Col ; /* pointer into A of Col [0..n_col] array */
1493  Int n_col2 ; /* number of non-dense, non-empty columns */
1494  Int n_row2 ; /* number of non-dense, non-empty rows */
1495  Int ngarbage ; /* number of garbage collections performed */
1496  Int max_deg ; /* maximum row degree */
1497  double default_knobs [COLAMD_KNOBS] ; /* default knobs array */
1498  Int aggressive ; /* do aggressive absorption */
1499  int ok ;
1500 
1501 #ifndef NDEBUG
1502  colamd_get_debug ("colamd") ;
1503 #endif /* NDEBUG */
1504 
1505  /* === Check the input arguments ======================================== */
1506 
1507  if (!stats)
1508  {
1509  DEBUG0 (("colamd: stats not present\n")) ;
1510  return (FALSE) ;
1511  }
1512  for (i = 0 ; i < COLAMD_STATS ; i++)
1513  {
1514  stats [i] = 0 ;
1515  }
1516  stats [COLAMD_STATUS] = COLAMD_OK ;
1517  stats [COLAMD_INFO1] = -1 ;
1518  stats [COLAMD_INFO2] = -1 ;
1519 
1520  if (!A) /* A is not present */
1521  {
1523  DEBUG0 (("colamd: A not present\n")) ;
1524  return (FALSE) ;
1525  }
1526 
1527  if (!p) /* p is not present */
1528  {
1530  DEBUG0 (("colamd: p not present\n")) ;
1531  return (FALSE) ;
1532  }
1533 
1534  if (n_row < 0) /* n_row must be >= 0 */
1535  {
1537  stats [COLAMD_INFO1] = n_row ;
1538  DEBUG0 (("colamd: nrow negative %d\n", n_row)) ;
1539  return (FALSE) ;
1540  }
1541 
1542  if (n_col < 0) /* n_col must be >= 0 */
1543  {
1545  stats [COLAMD_INFO1] = n_col ;
1546  DEBUG0 (("colamd: ncol negative %d\n", n_col)) ;
1547  return (FALSE) ;
1548  }
1549 
1550  nnz = p [n_col] ;
1551  if (nnz < 0) /* nnz must be >= 0 */
1552  {
1554  stats [COLAMD_INFO1] = nnz ;
1555  DEBUG0 (("colamd: number of entries negative %d\n", nnz)) ;
1556  return (FALSE) ;
1557  }
1558 
1559  if (p [0] != 0)
1560  {
1562  stats [COLAMD_INFO1] = p [0] ;
1563  DEBUG0 (("colamd: p[0] not zero %d\n", p [0])) ;
1564  return (FALSE) ;
1565  }
1566 
1567  /* === If no knobs, set default knobs =================================== */
1568 
1569  if (!knobs)
1570  {
1571  COLAMD_set_defaults (default_knobs) ;
1572  knobs = default_knobs ;
1573  }
1574 
1575  aggressive = (knobs [COLAMD_AGGRESSIVE] != FALSE) ;
1576 
1577  /* === Allocate the Row and Col arrays from array A ===================== */
1578 
1579  ok = TRUE ;
1580  Col_size = COLAMD_C (n_col, &ok) ; /* size of Col array of structs */
1581  Row_size = COLAMD_R (n_row, &ok) ; /* size of Row array of structs */
1582 
1583  /* need = 2*nnz + n_col + Col_size + Row_size ; */
1584  need = t_mult (nnz, 2, &ok) ;
1585  need = t_add (need, n_col, &ok) ;
1586  need = t_add (need, Col_size, &ok) ;
1587  need = t_add (need, Row_size, &ok) ;
1588 
1589  if (!ok || need > (size_t) Alen || need > Int_MAX)
1590  {
1591  /* not enough space in array A to perform the ordering */
1593  stats [COLAMD_INFO1] = need ;
1594  stats [COLAMD_INFO2] = Alen ;
1595  DEBUG0 (("colamd: Need Alen >= %d, given only Alen = %d\n", need,Alen));
1596  return (FALSE) ;
1597  }
1598 
1599  Alen -= Col_size + Row_size ;
1600  Col = (Colamd_Col *) &A [Alen] ;
1601  Row = (Colamd_Row *) &A [Alen + Col_size] ;
1602 
1603  /* === Construct the row and column data structures ===================== */
1604 
1605  if (!init_rows_cols (n_row, n_col, Row, Col, A, p, stats))
1606  {
1607  /* input matrix is invalid */
1608  DEBUG0 (("colamd: Matrix invalid\n")) ;
1609  return (FALSE) ;
1610  }
1611 
1612  /* === Initialize scores, kill dense rows/columns ======================= */
1613 
1614  init_scoring (n_row, n_col, Row, Col, A, p, knobs,
1615  &n_row2, &n_col2, &max_deg) ;
1616 
1617  /* === Order the supercolumns =========================================== */
1618 
1619  ngarbage = find_ordering (n_row, n_col, Alen, Row, Col, A, p,
1620  n_col2, max_deg, 2*nnz, aggressive) ;
1621 
1622  /* === Order the non-principal columns ================================== */
1623 
1624  order_children (n_col, Col, p) ;
1625 
1626  /* === Return statistics in stats ======================================= */
1627 
1628  stats [COLAMD_DENSE_ROW] = n_row - n_row2 ;
1629  stats [COLAMD_DENSE_COL] = n_col - n_col2 ;
1630  stats [COLAMD_DEFRAG_COUNT] = ngarbage ;
1631  DEBUG0 (("colamd: done.\n")) ;
1632  return (TRUE) ;
1633 }
1634 
1635 
1636 /* ========================================================================== */
1637 /* === colamd_report ======================================================== */
1638 /* ========================================================================== */
1639 
1640 PUBLIC void COLAMD_report
1642  Int stats [COLAMD_STATS]
1643 )
1644 {
1645  print_report ("colamd", stats) ;
1646 }
1647 
1648 
1649 /* ========================================================================== */
1650 /* === symamd_report ======================================================== */
1651 /* ========================================================================== */
1652 
1653 PUBLIC void SYMAMD_report
1655  Int stats [COLAMD_STATS]
1656 )
1657 {
1658  print_report ("symamd", stats) ;
1659 }
1660 
1661 
1662 
1663 /* ========================================================================== */
1664 /* === NON-USER-CALLABLE ROUTINES: ========================================== */
1665 /* ========================================================================== */
1666 
1667 /* There are no user-callable routines beyond this point in the file */
1668 
1669 
1670 /* ========================================================================== */
1671 /* === init_rows_cols ======================================================= */
1672 /* ========================================================================== */
1673 
1674 /*
1675  Takes the column form of the matrix in A and creates the row form of the
1676  matrix. Also, row and column attributes are stored in the Col and Row
1677  structs. If the columns are un-sorted or contain duplicate row indices,
1678  this routine will also sort and remove duplicate row indices from the
1679  column form of the matrix. Returns FALSE if the matrix is invalid,
1680  TRUE otherwise. Not user-callable.
1681 */
1682 
1683 PRIVATE Int init_rows_cols /* returns TRUE if OK, or FALSE otherwise */
1685  /* === Parameters ======================================================= */
1686 
1687  Int n_row, /* number of rows of A */
1688  Int n_col, /* number of columns of A */
1689  Colamd_Row Row [], /* of size n_row+1 */
1690  Colamd_Col Col [], /* of size n_col+1 */
1691  Int A [], /* row indices of A, of size Alen */
1692  Int p [], /* pointers to columns in A, of size n_col+1 */
1693  Int stats [COLAMD_STATS] /* colamd statistics */
1694 )
1695 {
1696  /* === Local variables ================================================== */
1697 
1698  Int col ; /* a column index */
1699  Int row ; /* a row index */
1700  Int *cp ; /* a column pointer */
1701  Int *cp_end ; /* a pointer to the end of a column */
1702  Int *rp ; /* a row pointer */
1703  Int *rp_end ; /* a pointer to the end of a row */
1704  Int last_row ; /* previous row */
1705 
1706  /* === Initialize columns, and check column pointers ==================== */
1707 
1708  for (col = 0 ; col < n_col ; col++)
1709  {
1710  Col [col].start = p [col] ;
1711  Col [col].length = p [col+1] - p [col] ;
1712 
1713  if (Col [col].length < 0)
1714  {
1715  /* column pointers must be non-decreasing */
1717  stats [COLAMD_INFO1] = col ;
1718  stats [COLAMD_INFO2] = Col [col].length ;
1719  DEBUG0 (("colamd: col %d length %d < 0\n", col, Col [col].length)) ;
1720  return (FALSE) ;
1721  }
1722 
1723  Col [col].shared1.thickness = 1 ;
1724  Col [col].shared2.score = 0 ;
1725  Col [col].shared3.prev = EMPTY ;
1726  Col [col].shared4.degree_next = EMPTY ;
1727  }
1728 
1729  /* p [0..n_col] no longer needed, used as "head" in subsequent routines */
1730 
1731  /* === Scan columns, compute row degrees, and check row indices ========= */
1732 
1733  stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
1734 
1735  for (row = 0 ; row < n_row ; row++)
1736  {
1737  Row [row].length = 0 ;
1738  Row [row].shared2.mark = -1 ;
1739  }
1740 
1741  for (col = 0 ; col < n_col ; col++)
1742  {
1743  last_row = -1 ;
1744 
1745  cp = &A [p [col]] ;
1746  cp_end = &A [p [col+1]] ;
1747 
1748  while (cp < cp_end)
1749  {
1750  row = *cp++ ;
1751 
1752  /* make sure row indices within range */
1753  if (row < 0 || row >= n_row)
1754  {
1756  stats [COLAMD_INFO1] = col ;
1757  stats [COLAMD_INFO2] = row ;
1758  stats [COLAMD_INFO3] = n_row ;
1759  DEBUG0 (("colamd: row %d col %d out of bounds\n", row, col)) ;
1760  return (FALSE) ;
1761  }
1762 
1763  if (row <= last_row || Row [row].shared2.mark == col)
1764  {
1765  /* row index are unsorted or repeated (or both), thus col */
1766  /* is jumbled. This is a notice, not an error condition. */
1768  stats [COLAMD_INFO1] = col ;
1769  stats [COLAMD_INFO2] = row ;
1770  (stats [COLAMD_INFO3]) ++ ;
1771  DEBUG1 (("colamd: row %d col %d unsorted/duplicate\n",row,col));
1772  }
1773 
1774  if (Row [row].shared2.mark != col)
1775  {
1776  Row [row].length++ ;
1777  }
1778  else
1779  {
1780  /* this is a repeated entry in the column, */
1781  /* it will be removed */
1782  Col [col].length-- ;
1783  }
1784 
1785  /* mark the row as having been seen in this column */
1786  Row [row].shared2.mark = col ;
1787 
1788  last_row = row ;
1789  }
1790  }
1791 
1792  /* === Compute row pointers ============================================= */
1793 
1794  /* row form of the matrix starts directly after the column */
1795  /* form of matrix in A */
1796  Row [0].start = p [n_col] ;
1797  Row [0].shared1.p = Row [0].start ;
1798  Row [0].shared2.mark = -1 ;
1799  for (row = 1 ; row < n_row ; row++)
1800  {
1801  Row [row].start = Row [row-1].start + Row [row-1].length ;
1802  Row [row].shared1.p = Row [row].start ;
1803  Row [row].shared2.mark = -1 ;
1804  }
1805 
1806  /* === Create row form ================================================== */
1807 
1808  if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1809  {
1810  /* if cols jumbled, watch for repeated row indices */
1811  for (col = 0 ; col < n_col ; col++)
1812  {
1813  cp = &A [p [col]] ;
1814  cp_end = &A [p [col+1]] ;
1815  while (cp < cp_end)
1816  {
1817  row = *cp++ ;
1818  if (Row [row].shared2.mark != col)
1819  {
1820  A [(Row [row].shared1.p)++] = col ;
1821  Row [row].shared2.mark = col ;
1822  }
1823  }
1824  }
1825  }
1826  else
1827  {
1828  /* if cols not jumbled, we don't need the mark (this is faster) */
1829  for (col = 0 ; col < n_col ; col++)
1830  {
1831  cp = &A [p [col]] ;
1832  cp_end = &A [p [col+1]] ;
1833  while (cp < cp_end)
1834  {
1835  A [(Row [*cp++].shared1.p)++] = col ;
1836  }
1837  }
1838  }
1839 
1840  /* === Clear the row marks and set row degrees ========================== */
1841 
1842  for (row = 0 ; row < n_row ; row++)
1843  {
1844  Row [row].shared2.mark = 0 ;
1845  Row [row].shared1.degree = Row [row].length ;
1846  }
1847 
1848  /* === See if we need to re-create columns ============================== */
1849 
1850  if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1851  {
1852  DEBUG0 (("colamd: reconstructing column form, matrix jumbled\n")) ;
1853 
1854 #ifndef NDEBUG
1855  /* make sure column lengths are correct */
1856  for (col = 0 ; col < n_col ; col++)
1857  {
1858  p [col] = Col [col].length ;
1859  }
1860  for (row = 0 ; row < n_row ; row++)
1861  {
1862  rp = &A [Row [row].start] ;
1863  rp_end = rp + Row [row].length ;
1864  while (rp < rp_end)
1865  {
1866  p [*rp++]-- ;
1867  }
1868  }
1869  for (col = 0 ; col < n_col ; col++)
1870  {
1871  ASSERT (p [col] == 0) ;
1872  }
1873  /* now p is all zero (different than when debugging is turned off) */
1874 #endif /* NDEBUG */
1875 
1876  /* === Compute col pointers ========================================= */
1877 
1878  /* col form of the matrix starts at A [0]. */
1879  /* Note, we may have a gap between the col form and the row */
1880  /* form if there were duplicate entries, if so, it will be */
1881  /* removed upon the first garbage collection */
1882  Col [0].start = 0 ;
1883  p [0] = Col [0].start ;
1884  for (col = 1 ; col < n_col ; col++)
1885  {
1886  /* note that the lengths here are for pruned columns, i.e. */
1887  /* no duplicate row indices will exist for these columns */
1888  Col [col].start = Col [col-1].start + Col [col-1].length ;
1889  p [col] = Col [col].start ;
1890  }
1891 
1892  /* === Re-create col form =========================================== */
1893 
1894  for (row = 0 ; row < n_row ; row++)
1895  {
1896  rp = &A [Row [row].start] ;
1897  rp_end = rp + Row [row].length ;
1898  while (rp < rp_end)
1899  {
1900  A [(p [*rp++])++] = row ;
1901  }
1902  }
1903  }
1904 
1905  /* === Done. Matrix is not (or no longer) jumbled ====================== */
1906 
1907  return (TRUE) ;
1908 }
1909 
1910 
1911 /* ========================================================================== */
1912 /* === init_scoring ========================================================= */
1913 /* ========================================================================== */
1914 
1915 /*
1916  Kills dense or empty columns and rows, calculates an initial score for
1917  each column, and places all columns in the degree lists. Not user-callable.
1918 */
1919 
1920 PRIVATE void init_scoring
1922  /* === Parameters ======================================================= */
1923 
1924  Int n_row, /* number of rows of A */
1925  Int n_col, /* number of columns of A */
1926  Colamd_Row Row [], /* of size n_row+1 */
1927  Colamd_Col Col [], /* of size n_col+1 */
1928  Int A [], /* column form and row form of A */
1929  Int head [], /* of size n_col+1 */
1930  double knobs [COLAMD_KNOBS],/* parameters */
1931  Int *p_n_row2, /* number of non-dense, non-empty rows */
1932  Int *p_n_col2, /* number of non-dense, non-empty columns */
1933  Int *p_max_deg /* maximum row degree */
1934 )
1935 {
1936  /* === Local variables ================================================== */
1937 
1938  Int c ; /* a column index */
1939  Int r, row ; /* a row index */
1940  Int *cp ; /* a column pointer */
1941  Int deg ; /* degree of a row or column */
1942  Int *cp_end ; /* a pointer to the end of a column */
1943  Int *new_cp ; /* new column pointer */
1944  Int col_length ; /* length of pruned column */
1945  Int score ; /* current column score */
1946  Int n_col2 ; /* number of non-dense, non-empty columns */
1947  Int n_row2 ; /* number of non-dense, non-empty rows */
1948  Int dense_row_count ; /* remove rows with more entries than this */
1949  Int dense_col_count ; /* remove cols with more entries than this */
1950  Int min_score ; /* smallest column score */
1951  Int max_deg ; /* maximum row degree */
1952  Int next_col ; /* Used to add to degree list.*/
1953 
1954 #ifndef NDEBUG
1955  Int debug_count ; /* debug only. */
1956 #endif /* NDEBUG */
1957 
1958  /* === Extract knobs ==================================================== */
1959 
1960  /* Note: if knobs contains a NaN, this is undefined: */
1961  if (knobs [COLAMD_DENSE_ROW] < 0)
1962  {
1963  /* only remove completely dense rows */
1964  dense_row_count = n_col-1 ;
1965  }
1966  else
1967  {
1968  dense_row_count = DENSE_DEGREE (knobs [COLAMD_DENSE_ROW], n_col) ;
1969  }
1970  if (knobs [COLAMD_DENSE_COL] < 0)
1971  {
1972  /* only remove completely dense columns */
1973  dense_col_count = n_row-1 ;
1974  }
1975  else
1976  {
1977  dense_col_count =
1978  DENSE_DEGREE (knobs [COLAMD_DENSE_COL], MIN (n_row, n_col)) ;
1979  }
1980 
1981  DEBUG1 (("colamd: densecount: %d %d\n", dense_row_count, dense_col_count)) ;
1982  max_deg = 0 ;
1983  n_col2 = n_col ;
1984  n_row2 = n_row ;
1985 
1986  /* === Kill empty columns =============================================== */
1987 
1988  /* Put the empty columns at the end in their natural order, so that LU */
1989  /* factorization can proceed as far as possible. */
1990  for (c = n_col-1 ; c >= 0 ; c--)
1991  {
1992  deg = Col [c].length ;
1993  if (deg == 0)
1994  {
1995  /* this is a empty column, kill and order it last */
1996  Col [c].shared2.order = --n_col2 ;
1997  KILL_PRINCIPAL_COL (c) ;
1998  }
1999  }
2000  DEBUG1 (("colamd: null columns killed: %d\n", n_col - n_col2)) ;
2001 
2002  /* === Kill dense columns =============================================== */
2003 
2004  /* Put the dense columns at the end, in their natural order */
2005  for (c = n_col-1 ; c >= 0 ; c--)
2006  {
2007  /* skip any dead columns */
2008  if (COL_IS_DEAD (c))
2009  {
2010  continue ;
2011  }
2012  deg = Col [c].length ;
2013  if (deg > dense_col_count)
2014  {
2015  /* this is a dense column, kill and order it last */
2016  Col [c].shared2.order = --n_col2 ;
2017  /* decrement the row degrees */
2018  cp = &A [Col [c].start] ;
2019  cp_end = cp + Col [c].length ;
2020  while (cp < cp_end)
2021  {
2022  Row [*cp++].shared1.degree-- ;
2023  }
2024  KILL_PRINCIPAL_COL (c) ;
2025  }
2026  }
2027  DEBUG1 (("colamd: Dense and null columns killed: %d\n", n_col - n_col2)) ;
2028 
2029  /* === Kill dense and empty rows ======================================== */
2030 
2031  for (r = 0 ; r < n_row ; r++)
2032  {
2033  deg = Row [r].shared1.degree ;
2034  ASSERT (deg >= 0 && deg <= n_col) ;
2035  if (deg > dense_row_count || deg == 0)
2036  {
2037  /* kill a dense or empty row */
2038  KILL_ROW (r) ;
2039  --n_row2 ;
2040  }
2041  else
2042  {
2043  /* keep track of max degree of remaining rows */
2044  max_deg = MAX (max_deg, deg) ;
2045  }
2046  }
2047  DEBUG1 (("colamd: Dense and null rows killed: %d\n", n_row - n_row2)) ;
2048 
2049  /* === Compute initial column scores ==================================== */
2050 
2051  /* At this point the row degrees are accurate. They reflect the number */
2052  /* of "live" (non-dense) columns in each row. No empty rows exist. */
2053  /* Some "live" columns may contain only dead rows, however. These are */
2054  /* pruned in the code below. */
2055 
2056  /* now find the initial matlab score for each column */
2057  for (c = n_col-1 ; c >= 0 ; c--)
2058  {
2059  /* skip dead column */
2060  if (COL_IS_DEAD (c))
2061  {
2062  continue ;
2063  }
2064  score = 0 ;
2065  cp = &A [Col [c].start] ;
2066  new_cp = cp ;
2067  cp_end = cp + Col [c].length ;
2068  while (cp < cp_end)
2069  {
2070  /* get a row */
2071  row = *cp++ ;
2072  /* skip if dead */
2073  if (ROW_IS_DEAD (row))
2074  {
2075  continue ;
2076  }
2077  /* compact the column */
2078  *new_cp++ = row ;
2079  /* add row's external degree */
2080  score += Row [row].shared1.degree - 1 ;
2081  /* guard against integer overflow */
2082  score = MIN (score, n_col) ;
2083  }
2084  /* determine pruned column length */
2085  col_length = (Int) (new_cp - &A [Col [c].start]) ;
2086  if (col_length == 0)
2087  {
2088  /* a newly-made null column (all rows in this col are "dense" */
2089  /* and have already been killed) */
2090  DEBUG2 (("Newly null killed: %d\n", c)) ;
2091  Col [c].shared2.order = --n_col2 ;
2092  KILL_PRINCIPAL_COL (c) ;
2093  }
2094  else
2095  {
2096  /* set column length and set score */
2097  ASSERT (score >= 0) ;
2098  ASSERT (score <= n_col) ;
2099  Col [c].length = col_length ;
2100  Col [c].shared2.score = score ;
2101  }
2102  }
2103  DEBUG1 (("colamd: Dense, null, and newly-null columns killed: %d\n",
2104  n_col-n_col2)) ;
2105 
2106  /* At this point, all empty rows and columns are dead. All live columns */
2107  /* are "clean" (containing no dead rows) and simplicial (no supercolumns */
2108  /* yet). Rows may contain dead columns, but all live rows contain at */
2109  /* least one live column. */
2110 
2111 #ifndef NDEBUG
2112  debug_structures (n_row, n_col, Row, Col, A, n_col2) ;
2113 #endif /* NDEBUG */
2114 
2115  /* === Initialize degree lists ========================================== */
2116 
2117 #ifndef NDEBUG
2118  debug_count = 0 ;
2119 #endif /* NDEBUG */
2120 
2121  /* clear the hash buckets */
2122  for (c = 0 ; c <= n_col ; c++)
2123  {
2124  head [c] = EMPTY ;
2125  }
2126  min_score = n_col ;
2127  /* place in reverse order, so low column indices are at the front */
2128  /* of the lists. This is to encourage natural tie-breaking */
2129  for (c = n_col-1 ; c >= 0 ; c--)
2130  {
2131  /* only add principal columns to degree lists */
2132  if (COL_IS_ALIVE (c))
2133  {
2134  DEBUG4 (("place %d score %d minscore %d ncol %d\n",
2135  c, Col [c].shared2.score, min_score, n_col)) ;
2136 
2137  /* === Add columns score to DList =============================== */
2138 
2139  score = Col [c].shared2.score ;
2140 
2141  ASSERT (min_score >= 0) ;
2142  ASSERT (min_score <= n_col) ;
2143  ASSERT (score >= 0) ;
2144  ASSERT (score <= n_col) ;
2145  ASSERT (head [score] >= EMPTY) ;
2146 
2147  /* now add this column to dList at proper score location */
2148  next_col = head [score] ;
2149  Col [c].shared3.prev = EMPTY ;
2150  Col [c].shared4.degree_next = next_col ;
2151 
2152  /* if there already was a column with the same score, set its */
2153  /* previous pointer to this new column */
2154  if (next_col != EMPTY)
2155  {
2156  Col [next_col].shared3.prev = c ;
2157  }
2158  head [score] = c ;
2159 
2160  /* see if this score is less than current min */
2161  min_score = MIN (min_score, score) ;
2162 
2163 #ifndef NDEBUG
2164  debug_count++ ;
2165 #endif /* NDEBUG */
2166 
2167  }
2168  }
2169 
2170 #ifndef NDEBUG
2171  DEBUG1 (("colamd: Live cols %d out of %d, non-princ: %d\n",
2172  debug_count, n_col, n_col-debug_count)) ;
2173  ASSERT (debug_count == n_col2) ;
2174  debug_deg_lists (n_row, n_col, Row, Col, head, min_score, n_col2, max_deg) ;
2175 #endif /* NDEBUG */
2176 
2177  /* === Return number of remaining columns, and max row degree =========== */
2178 
2179  *p_n_col2 = n_col2 ;
2180  *p_n_row2 = n_row2 ;
2181  *p_max_deg = max_deg ;
2182 }
2183 
2184 
2185 /* ========================================================================== */
2186 /* === find_ordering ======================================================== */
2187 /* ========================================================================== */
2188 
2189 /*
2190  Order the principal columns of the supercolumn form of the matrix
2191  (no supercolumns on input). Uses a minimum approximate column minimum
2192  degree ordering method. Not user-callable.
2193 */
2194 
2195 PRIVATE Int find_ordering /* return the number of garbage collections */
2197  /* === Parameters ======================================================= */
2198 
2199  Int n_row, /* number of rows of A */
2200  Int n_col, /* number of columns of A */
2201  Int Alen, /* size of A, 2*nnz + n_col or larger */
2202  Colamd_Row Row [], /* of size n_row+1 */
2203  Colamd_Col Col [], /* of size n_col+1 */
2204  Int A [], /* column form and row form of A */
2205  Int head [], /* of size n_col+1 */
2206  Int n_col2, /* Remaining columns to order */
2207  Int max_deg, /* Maximum row degree */
2208  Int pfree, /* index of first free slot (2*nnz on entry) */
2209  Int aggressive
2210 )
2211 {
2212  /* === Local variables ================================================== */
2213 
2214  Int k ; /* current pivot ordering step */
2215  Int pivot_col ; /* current pivot column */
2216  Int *cp ; /* a column pointer */
2217  Int *rp ; /* a row pointer */
2218  Int pivot_row ; /* current pivot row */
2219  Int *new_cp ; /* modified column pointer */
2220  Int *new_rp ; /* modified row pointer */
2221  Int pivot_row_start ; /* pointer to start of pivot row */
2222  Int pivot_row_degree ; /* number of columns in pivot row */
2223  Int pivot_row_length ; /* number of supercolumns in pivot row */
2224  Int pivot_col_score ; /* score of pivot column */
2225  Int needed_memory ; /* free space needed for pivot row */
2226  Int *cp_end ; /* pointer to the end of a column */
2227  Int *rp_end ; /* pointer to the end of a row */
2228  Int row ; /* a row index */
2229  Int col ; /* a column index */
2230  Int max_score ; /* maximum possible score */
2231  Int cur_score ; /* score of current column */
2232  unsigned Int hash ; /* hash value for supernode detection */
2233  Int head_column ; /* head of hash bucket */
2234  Int first_col ; /* first column in hash bucket */
2235  Int tag_mark ; /* marker value for mark array */
2236  Int row_mark ; /* Row [row].shared2.mark */
2237  Int set_difference ; /* set difference size of row with pivot row */
2238  Int min_score ; /* smallest column score */
2239  Int col_thickness ; /* "thickness" (no. of columns in a supercol) */
2240  Int max_mark ; /* maximum value of tag_mark */
2241  Int pivot_col_thickness ; /* number of columns represented by pivot col */
2242  Int prev_col ; /* Used by Dlist operations. */
2243  Int next_col ; /* Used by Dlist operations. */
2244  Int ngarbage ; /* number of garbage collections performed */
2245 
2246 #ifndef NDEBUG
2247  Int debug_d ; /* debug loop counter */
2248  Int debug_step = 0 ; /* debug loop counter */
2249 #endif /* NDEBUG */
2250 
2251  /* === Initialization and clear mark ==================================== */
2252 
2253  max_mark = INT_MAX - n_col ; /* INT_MAX defined in <limits.h> */
2254  tag_mark = clear_mark (0, max_mark, n_row, Row) ;
2255  min_score = 0 ;
2256  ngarbage = 0 ;
2257  DEBUG1 (("colamd: Ordering, n_col2=%d\n", n_col2)) ;
2258 
2259  /* === Order the columns ================================================ */
2260 
2261  for (k = 0 ; k < n_col2 ; /* 'k' is incremented below */)
2262  {
2263 
2264 #ifndef NDEBUG
2265  if (debug_step % 100 == 0)
2266  {
2267  DEBUG2 (("\n... Step k: %d out of n_col2: %d\n", k, n_col2)) ;
2268  }
2269  else
2270  {
2271  DEBUG3 (("\n----------Step k: %d out of n_col2: %d\n", k, n_col2)) ;
2272  }
2273  debug_step++ ;
2274  debug_deg_lists (n_row, n_col, Row, Col, head,
2275  min_score, n_col2-k, max_deg) ;
2276  debug_matrix (n_row, n_col, Row, Col, A) ;
2277 #endif /* NDEBUG */
2278 
2279  /* === Select pivot column, and order it ============================ */
2280 
2281  /* make sure degree list isn't empty */
2282  ASSERT (min_score >= 0) ;
2283  ASSERT (min_score <= n_col) ;
2284  ASSERT (head [min_score] >= EMPTY) ;
2285 
2286 #ifndef NDEBUG
2287  for (debug_d = 0 ; debug_d < min_score ; debug_d++)
2288  {
2289  ASSERT (head [debug_d] == EMPTY) ;
2290  }
2291 #endif /* NDEBUG */
2292 
2293  /* get pivot column from head of minimum degree list */
2294  while (head [min_score] == EMPTY && min_score < n_col)
2295  {
2296  min_score++ ;
2297  }
2298  pivot_col = head [min_score] ;
2299  ASSERT (pivot_col >= 0 && pivot_col <= n_col) ;
2300  next_col = Col [pivot_col].shared4.degree_next ;
2301  head [min_score] = next_col ;
2302  if (next_col != EMPTY)
2303  {
2304  Col [next_col].shared3.prev = EMPTY ;
2305  }
2306 
2307  ASSERT (COL_IS_ALIVE (pivot_col)) ;
2308 
2309  /* remember score for defrag check */
2310  pivot_col_score = Col [pivot_col].shared2.score ;
2311 
2312  /* the pivot column is the kth column in the pivot order */
2313  Col [pivot_col].shared2.order = k ;
2314 
2315  /* increment order count by column thickness */
2316  pivot_col_thickness = Col [pivot_col].shared1.thickness ;
2317  k += pivot_col_thickness ;
2318  ASSERT (pivot_col_thickness > 0) ;
2319  DEBUG3 (("Pivot col: %d thick %d\n", pivot_col, pivot_col_thickness)) ;
2320 
2321  /* === Garbage_collection, if necessary ============================= */
2322 
2323  needed_memory = MIN (pivot_col_score, n_col - k) ;
2324  if (pfree + needed_memory >= Alen)
2325  {
2326  pfree = garbage_collection (n_row, n_col, Row, Col, A, &A [pfree]) ;
2327  ngarbage++ ;
2328  /* after garbage collection we will have enough */
2329  ASSERT (pfree + needed_memory < Alen) ;
2330  /* garbage collection has wiped out the Row[].shared2.mark array */
2331  tag_mark = clear_mark (0, max_mark, n_row, Row) ;
2332 
2333 #ifndef NDEBUG
2334  debug_matrix (n_row, n_col, Row, Col, A) ;
2335 #endif /* NDEBUG */
2336  }
2337 
2338  /* === Compute pivot row pattern ==================================== */
2339 
2340  /* get starting location for this new merged row */
2341  pivot_row_start = pfree ;
2342 
2343  /* initialize new row counts to zero */
2344  pivot_row_degree = 0 ;
2345 
2346  /* tag pivot column as having been visited so it isn't included */
2347  /* in merged pivot row */
2348  Col [pivot_col].shared1.thickness = -pivot_col_thickness ;
2349 
2350  /* pivot row is the union of all rows in the pivot column pattern */
2351  cp = &A [Col [pivot_col].start] ;
2352  cp_end = cp + Col [pivot_col].length ;
2353  while (cp < cp_end)
2354  {
2355  /* get a row */
2356  row = *cp++ ;
2357  DEBUG4 (("Pivot col pattern %d %d\n", ROW_IS_ALIVE (row), row)) ;
2358  /* skip if row is dead */
2359  if (ROW_IS_ALIVE (row))
2360  {
2361  rp = &A [Row [row].start] ;
2362  rp_end = rp + Row [row].length ;
2363  while (rp < rp_end)
2364  {
2365  /* get a column */
2366  col = *rp++ ;
2367  /* add the column, if alive and untagged */
2368  col_thickness = Col [col].shared1.thickness ;
2369  if (col_thickness > 0 && COL_IS_ALIVE (col))
2370  {
2371  /* tag column in pivot row */
2372  Col [col].shared1.thickness = -col_thickness ;
2373  ASSERT (pfree < Alen) ;
2374  /* place column in pivot row */
2375  A [pfree++] = col ;
2376  pivot_row_degree += col_thickness ;
2377  }
2378  }
2379  }
2380  }
2381 
2382  /* clear tag on pivot column */
2383  Col [pivot_col].shared1.thickness = pivot_col_thickness ;
2384  max_deg = MAX (max_deg, pivot_row_degree) ;
2385 
2386 #ifndef NDEBUG
2387  DEBUG3 (("check2\n")) ;
2388  debug_mark (n_row, Row, tag_mark, max_mark) ;
2389 #endif /* NDEBUG */
2390 
2391  /* === Kill all rows used to construct pivot row ==================== */
2392 
2393  /* also kill pivot row, temporarily */
2394  cp = &A [Col [pivot_col].start] ;
2395  cp_end = cp + Col [pivot_col].length ;
2396  while (cp < cp_end)
2397  {
2398  /* may be killing an already dead row */
2399  row = *cp++ ;
2400  DEBUG3 (("Kill row in pivot col: %d\n", row)) ;
2401  KILL_ROW (row) ;
2402  }
2403 
2404  /* === Select a row index to use as the new pivot row =============== */
2405 
2406  pivot_row_length = pfree - pivot_row_start ;
2407  if (pivot_row_length > 0)
2408  {
2409  /* pick the "pivot" row arbitrarily (first row in col) */
2410  pivot_row = A [Col [pivot_col].start] ;
2411  DEBUG3 (("Pivotal row is %d\n", pivot_row)) ;
2412  }
2413  else
2414  {
2415  /* there is no pivot row, since it is of zero length */
2416  pivot_row = EMPTY ;
2417  ASSERT (pivot_row_length == 0) ;
2418  }
2419  ASSERT (Col [pivot_col].length > 0 || pivot_row_length == 0) ;
2420 
2421  /* === Approximate degree computation =============================== */
2422 
2423  /* Here begins the computation of the approximate degree. The column */
2424  /* score is the sum of the pivot row "length", plus the size of the */
2425  /* set differences of each row in the column minus the pattern of the */
2426  /* pivot row itself. The column ("thickness") itself is also */
2427  /* excluded from the column score (we thus use an approximate */
2428  /* external degree). */
2429 
2430  /* The time taken by the following code (compute set differences, and */
2431  /* add them up) is proportional to the size of the data structure */
2432  /* being scanned - that is, the sum of the sizes of each column in */
2433  /* the pivot row. Thus, the amortized time to compute a column score */
2434  /* is proportional to the size of that column (where size, in this */
2435  /* context, is the column "length", or the number of row indices */
2436  /* in that column). The number of row indices in a column is */
2437  /* monotonically non-decreasing, from the length of the original */
2438  /* column on input to colamd. */
2439 
2440  /* === Compute set differences ====================================== */
2441 
2442  DEBUG3 (("** Computing set differences phase. **\n")) ;
2443 
2444  /* pivot row is currently dead - it will be revived later. */
2445 
2446  DEBUG3 (("Pivot row: ")) ;
2447  /* for each column in pivot row */
2448  rp = &A [pivot_row_start] ;
2449  rp_end = rp + pivot_row_length ;
2450  while (rp < rp_end)
2451  {
2452  col = *rp++ ;
2453  ASSERT (COL_IS_ALIVE (col) && col != pivot_col) ;
2454  DEBUG3 (("Col: %d\n", col)) ;
2455 
2456  /* clear tags used to construct pivot row pattern */
2457  col_thickness = -Col [col].shared1.thickness ;
2458  ASSERT (col_thickness > 0) ;
2459  Col [col].shared1.thickness = col_thickness ;
2460 
2461  /* === Remove column from degree list =========================== */
2462 
2463  cur_score = Col [col].shared2.score ;
2464  prev_col = Col [col].shared3.prev ;
2465  next_col = Col [col].shared4.degree_next ;
2466  ASSERT (cur_score >= 0) ;
2467  ASSERT (cur_score <= n_col) ;
2468  ASSERT (cur_score >= EMPTY) ;
2469  if (prev_col == EMPTY)
2470  {
2471  head [cur_score] = next_col ;
2472  }
2473  else
2474  {
2475  Col [prev_col].shared4.degree_next = next_col ;
2476  }
2477  if (next_col != EMPTY)
2478  {
2479  Col [next_col].shared3.prev = prev_col ;
2480  }
2481 
2482  /* === Scan the column ========================================== */
2483 
2484  cp = &A [Col [col].start] ;
2485  cp_end = cp + Col [col].length ;
2486  while (cp < cp_end)
2487  {
2488  /* get a row */
2489  row = *cp++ ;
2490  row_mark = Row [row].shared2.mark ;
2491  /* skip if dead */
2492  if (ROW_IS_MARKED_DEAD (row_mark))
2493  {
2494  continue ;
2495  }
2496  ASSERT (row != pivot_row) ;
2497  set_difference = row_mark - tag_mark ;
2498  /* check if the row has been seen yet */
2499  if (set_difference < 0)
2500  {
2501  ASSERT (Row [row].shared1.degree <= max_deg) ;
2502  set_difference = Row [row].shared1.degree ;
2503  }
2504  /* subtract column thickness from this row's set difference */
2505  set_difference -= col_thickness ;
2506  ASSERT (set_difference >= 0) ;
2507  /* absorb this row if the set difference becomes zero */
2508  if (set_difference == 0 && aggressive)
2509  {
2510  DEBUG3 (("aggressive absorption. Row: %d\n", row)) ;
2511  KILL_ROW (row) ;
2512  }
2513  else
2514  {
2515  /* save the new mark */
2516  Row [row].shared2.mark = set_difference + tag_mark ;
2517  }
2518  }
2519  }
2520 
2521 #ifndef NDEBUG
2522  debug_deg_lists (n_row, n_col, Row, Col, head,
2523  min_score, n_col2-k-pivot_row_degree, max_deg) ;
2524 #endif /* NDEBUG */
2525 
2526  /* === Add up set differences for each column ======================= */
2527 
2528  DEBUG3 (("** Adding set differences phase. **\n")) ;
2529 
2530  /* for each column in pivot row */
2531  rp = &A [pivot_row_start] ;
2532  rp_end = rp + pivot_row_length ;
2533  while (rp < rp_end)
2534  {
2535  /* get a column */
2536  col = *rp++ ;
2537  ASSERT (COL_IS_ALIVE (col) && col != pivot_col) ;
2538  hash = 0 ;
2539  cur_score = 0 ;
2540  cp = &A [Col [col].start] ;
2541  /* compact the column */
2542  new_cp = cp ;
2543  cp_end = cp + Col [col].length ;
2544 
2545  DEBUG4 (("Adding set diffs for Col: %d.\n", col)) ;
2546 
2547  while (cp < cp_end)
2548  {
2549  /* get a row */
2550  row = *cp++ ;
2551  ASSERT(row >= 0 && row < n_row) ;
2552  row_mark = Row [row].shared2.mark ;
2553  /* skip if dead */
2554  if (ROW_IS_MARKED_DEAD (row_mark))
2555  {
2556  DEBUG4 ((" Row %d, dead\n", row)) ;
2557  continue ;
2558  }
2559  DEBUG4 ((" Row %d, set diff %d\n", row, row_mark-tag_mark));
2560  ASSERT (row_mark >= tag_mark) ;
2561  /* compact the column */
2562  *new_cp++ = row ;
2563  /* compute hash function */
2564  hash += row ;
2565  /* add set difference */
2566  cur_score += row_mark - tag_mark ;
2567  /* integer overflow... */
2568  cur_score = MIN (cur_score, n_col) ;
2569  }
2570 
2571  /* recompute the column's length */
2572  Col [col].length = (Int) (new_cp - &A [Col [col].start]) ;
2573 
2574  /* === Further mass elimination ================================= */
2575 
2576  if (Col [col].length == 0)
2577  {
2578  DEBUG4 (("further mass elimination. Col: %d\n", col)) ;
2579  /* nothing left but the pivot row in this column */
2580  KILL_PRINCIPAL_COL (col) ;
2581  pivot_row_degree -= Col [col].shared1.thickness ;
2582  ASSERT (pivot_row_degree >= 0) ;
2583  /* order it */
2584  Col [col].shared2.order = k ;
2585  /* increment order count by column thickness */
2586  k += Col [col].shared1.thickness ;
2587  }
2588  else
2589  {
2590  /* === Prepare for supercolumn detection ==================== */
2591 
2592  DEBUG4 (("Preparing supercol detection for Col: %d.\n", col)) ;
2593 
2594  /* save score so far */
2595  Col [col].shared2.score = cur_score ;
2596 
2597  /* add column to hash table, for supercolumn detection */
2598  hash %= n_col + 1 ;
2599 
2600  DEBUG4 ((" Hash = %d, n_col = %d.\n", hash, n_col)) ;
2601  ASSERT (((Int) hash) <= n_col) ;
2602 
2603  head_column = head [hash] ;
2604  if (head_column > EMPTY)
2605  {
2606  /* degree list "hash" is non-empty, use prev (shared3) of */
2607  /* first column in degree list as head of hash bucket */
2608  first_col = Col [head_column].shared3.headhash ;
2609  Col [head_column].shared3.headhash = col ;
2610  }
2611  else
2612  {
2613  /* degree list "hash" is empty, use head as hash bucket */
2614  first_col = - (head_column + 2) ;
2615  head [hash] = - (col + 2) ;
2616  }
2617  Col [col].shared4.hash_next = first_col ;
2618 
2619  /* save hash function in Col [col].shared3.hash */
2620  Col [col].shared3.hash = (Int) hash ;
2621  ASSERT (COL_IS_ALIVE (col)) ;
2622  }
2623  }
2624 
2625  /* The approximate external column degree is now computed. */
2626 
2627  /* === Supercolumn detection ======================================== */
2628 
2629  DEBUG3 (("** Supercolumn detection phase. **\n")) ;
2630 
2632 
2633 #ifndef NDEBUG
2634  n_col, Row,
2635 #endif /* NDEBUG */
2636 
2637  Col, A, head, pivot_row_start, pivot_row_length) ;
2638 
2639  /* === Kill the pivotal column ====================================== */
2640 
2641  KILL_PRINCIPAL_COL (pivot_col) ;
2642 
2643  /* === Clear mark =================================================== */
2644 
2645  tag_mark = clear_mark (tag_mark+max_deg+1, max_mark, n_row, Row) ;
2646 
2647 #ifndef NDEBUG
2648  DEBUG3 (("check3\n")) ;
2649  debug_mark (n_row, Row, tag_mark, max_mark) ;
2650 #endif /* NDEBUG */
2651 
2652  /* === Finalize the new pivot row, and column scores ================ */
2653 
2654  DEBUG3 (("** Finalize scores phase. **\n")) ;
2655 
2656  /* for each column in pivot row */
2657  rp = &A [pivot_row_start] ;
2658  /* compact the pivot row */
2659  new_rp = rp ;
2660  rp_end = rp + pivot_row_length ;
2661  while (rp < rp_end)
2662  {
2663  col = *rp++ ;
2664  /* skip dead columns */
2665  if (COL_IS_DEAD (col))
2666  {
2667  continue ;
2668  }
2669  *new_rp++ = col ;
2670  /* add new pivot row to column */
2671  A [Col [col].start + (Col [col].length++)] = pivot_row ;
2672 
2673  /* retrieve score so far and add on pivot row's degree. */
2674  /* (we wait until here for this in case the pivot */
2675  /* row's degree was reduced due to mass elimination). */
2676  cur_score = Col [col].shared2.score + pivot_row_degree ;
2677 
2678  /* calculate the max possible score as the number of */
2679  /* external columns minus the 'k' value minus the */
2680  /* columns thickness */
2681  max_score = n_col - k - Col [col].shared1.thickness ;
2682 
2683  /* make the score the external degree of the union-of-rows */
2684  cur_score -= Col [col].shared1.thickness ;
2685 
2686  /* make sure score is less or equal than the max score */
2687  cur_score = MIN (cur_score, max_score) ;
2688  ASSERT (cur_score >= 0) ;
2689 
2690  /* store updated score */
2691  Col [col].shared2.score = cur_score ;
2692 
2693  /* === Place column back in degree list ========================= */
2694 
2695  ASSERT (min_score >= 0) ;
2696  ASSERT (min_score <= n_col) ;
2697  ASSERT (cur_score >= 0) ;
2698  ASSERT (cur_score <= n_col) ;
2699  ASSERT (head [cur_score] >= EMPTY) ;
2700  next_col = head [cur_score] ;
2701  Col [col].shared4.degree_next = next_col ;
2702  Col [col].shared3.prev = EMPTY ;
2703  if (next_col != EMPTY)
2704  {
2705  Col [next_col].shared3.prev = col ;
2706  }
2707  head [cur_score] = col ;
2708 
2709  /* see if this score is less than current min */
2710  min_score = MIN (min_score, cur_score) ;
2711 
2712  }
2713 
2714 #ifndef NDEBUG
2715  debug_deg_lists (n_row, n_col, Row, Col, head,
2716  min_score, n_col2-k, max_deg) ;
2717 #endif /* NDEBUG */
2718 
2719  /* === Resurrect the new pivot row ================================== */
2720 
2721  if (pivot_row_degree > 0)
2722  {
2723  /* update pivot row length to reflect any cols that were killed */
2724  /* during super-col detection and mass elimination */
2725  Row [pivot_row].start = pivot_row_start ;
2726  Row [pivot_row].length = (Int) (new_rp - &A[pivot_row_start]) ;
2727  ASSERT (Row [pivot_row].length > 0) ;
2728  Row [pivot_row].shared1.degree = pivot_row_degree ;
2729  Row [pivot_row].shared2.mark = 0 ;
2730  /* pivot row is no longer dead */
2731 
2732  DEBUG1 (("Resurrect Pivot_row %d deg: %d\n",
2733  pivot_row, pivot_row_degree)) ;
2734  }
2735  }
2736 
2737  /* === All principal columns have now been ordered ====================== */
2738 
2739  return (ngarbage) ;
2740 }
2741 
2742 
2743 /* ========================================================================== */
2744 /* === order_children ======================================================= */
2745 /* ========================================================================== */
2746 
2747 /*
2748  The find_ordering routine has ordered all of the principal columns (the
2749  representatives of the supercolumns). The non-principal columns have not
2750  yet been ordered. This routine orders those columns by walking up the
2751  parent tree (a column is a child of the column which absorbed it). The
2752  final permutation vector is then placed in p [0 ... n_col-1], with p [0]
2753  being the first column, and p [n_col-1] being the last. It doesn't look
2754  like it at first glance, but be assured that this routine takes time linear
2755  in the number of columns. Although not immediately obvious, the time
2756  taken by this routine is O (n_col), that is, linear in the number of
2757  columns. Not user-callable.
2758 */
2759 
2762  /* === Parameters ======================================================= */
2763 
2764  Int n_col, /* number of columns of A */
2765  Colamd_Col Col [], /* of size n_col+1 */
2766  Int p [] /* p [0 ... n_col-1] is the column permutation*/
2767 )
2768 {
2769  /* === Local variables ================================================== */
2770 
2771  Int i ; /* loop counter for all columns */
2772  Int c ; /* column index */
2773  Int parent ; /* index of column's parent */
2774  Int order ; /* column's order */
2775 
2776  /* === Order each non-principal column ================================== */
2777 
2778  for (i = 0 ; i < n_col ; i++)
2779  {
2780  /* find an un-ordered non-principal column */
2781  ASSERT (COL_IS_DEAD (i)) ;
2782  if (!COL_IS_DEAD_PRINCIPAL (i) && Col [i].shared2.order == EMPTY)
2783  {
2784  parent = i ;
2785  /* once found, find its principal parent */
2786  do
2787  {
2788  parent = Col [parent].shared1.parent ;
2789  } while (!COL_IS_DEAD_PRINCIPAL (parent)) ;
2790 
2791  /* now, order all un-ordered non-principal columns along path */
2792  /* to this parent. collapse tree at the same time */
2793  c = i ;
2794  /* get order of parent */
2795  order = Col [parent].shared2.order ;
2796 
2797  do
2798  {
2799  ASSERT (Col [c].shared2.order == EMPTY) ;
2800 
2801  /* order this column */
2802  Col [c].shared2.order = order++ ;
2803  /* collaps tree */
2804  Col [c].shared1.parent = parent ;
2805 
2806  /* get immediate parent of this column */
2807  c = Col [c].shared1.parent ;
2808 
2809  /* continue until we hit an ordered column. There are */
2810  /* guarranteed not to be anymore unordered columns */
2811  /* above an ordered column */
2812  } while (Col [c].shared2.order == EMPTY) ;
2813 
2814  /* re-order the super_col parent to largest order for this group */
2815  Col [parent].shared2.order = order ;
2816  }
2817  }
2818 
2819  /* === Generate the permutation ========================================= */
2820 
2821  for (c = 0 ; c < n_col ; c++)
2822  {
2823  p [Col [c].shared2.order] = c ;
2824  }
2825 }
2826 
2827 
2828 /* ========================================================================== */
2829 /* === detect_super_cols ==================================================== */
2830 /* ========================================================================== */
2831 
2832 /*
2833  Detects supercolumns by finding matches between columns in the hash buckets.
2834  Check amongst columns in the set A [row_start ... row_start + row_length-1].
2835  The columns under consideration are currently *not* in the degree lists,
2836  and have already been placed in the hash buckets.
2837 
2838  The hash bucket for columns whose hash function is equal to h is stored
2839  as follows:
2840 
2841  if head [h] is >= 0, then head [h] contains a degree list, so:
2842 
2843  head [h] is the first column in degree bucket h.
2844  Col [head [h]].headhash gives the first column in hash bucket h.
2845 
2846  otherwise, the degree list is empty, and:
2847 
2848  -(head [h] + 2) is the first column in hash bucket h.
2849 
2850  For a column c in a hash bucket, Col [c].shared3.prev is NOT a "previous
2851  column" pointer. Col [c].shared3.hash is used instead as the hash number
2852  for that column. The value of Col [c].shared4.hash_next is the next column
2853  in the same hash bucket.
2854 
2855  Assuming no, or "few" hash collisions, the time taken by this routine is
2856  linear in the sum of the sizes (lengths) of each column whose score has
2857  just been computed in the approximate degree computation.
2858  Not user-callable.
2859 */
2860 
2863  /* === Parameters ======================================================= */
2864 
2865 #ifndef NDEBUG
2866  /* these two parameters are only needed when debugging is enabled: */
2867  Int n_col, /* number of columns of A */
2868  Colamd_Row Row [], /* of size n_row+1 */
2869 #endif /* NDEBUG */
2870 
2871  Colamd_Col Col [], /* of size n_col+1 */
2872  Int A [], /* row indices of A */
2873  Int head [], /* head of degree lists and hash buckets */
2874  Int row_start, /* pointer to set of columns to check */
2875  Int row_length /* number of columns to check */
2876 )
2877 {
2878  /* === Local variables ================================================== */
2879 
2880  Int hash ; /* hash value for a column */
2881  Int *rp ; /* pointer to a row */
2882  Int c ; /* a column index */
2883  Int super_c ; /* column index of the column to absorb into */
2884  Int *cp1 ; /* column pointer for column super_c */
2885  Int *cp2 ; /* column pointer for column c */
2886  Int length ; /* length of column super_c */
2887  Int prev_c ; /* column preceding c in hash bucket */
2888  Int i ; /* loop counter */
2889  Int *rp_end ; /* pointer to the end of the row */
2890  Int col ; /* a column index in the row to check */
2891  Int head_column ; /* first column in hash bucket or degree list */
2892  Int first_col ; /* first column in hash bucket */
2893 
2894  /* === Consider each column in the row ================================== */
2895 
2896  rp = &A [row_start] ;
2897  rp_end = rp + row_length ;
2898  while (rp < rp_end)
2899  {
2900  col = *rp++ ;
2901  if (COL_IS_DEAD (col))
2902  {
2903  continue ;
2904  }
2905 
2906  /* get hash number for this column */
2907  hash = Col [col].shared3.hash ;
2908  ASSERT (hash <= n_col) ;
2909 
2910  /* === Get the first column in this hash bucket ===================== */
2911 
2912  head_column = head [hash] ;
2913  if (head_column > EMPTY)
2914  {
2915  first_col = Col [head_column].shared3.headhash ;
2916  }
2917  else
2918  {
2919  first_col = - (head_column + 2) ;
2920  }
2921 
2922  /* === Consider each column in the hash bucket ====================== */
2923 
2924  for (super_c = first_col ; super_c != EMPTY ;
2925  super_c = Col [super_c].shared4.hash_next)
2926  {
2927  ASSERT (COL_IS_ALIVE (super_c)) ;
2928  ASSERT (Col [super_c].shared3.hash == hash) ;
2929  length = Col [super_c].length ;
2930 
2931  /* prev_c is the column preceding column c in the hash bucket */
2932  prev_c = super_c ;
2933 
2934  /* === Compare super_c with all columns after it ================ */
2935 
2936  for (c = Col [super_c].shared4.hash_next ;
2937  c != EMPTY ; c = Col [c].shared4.hash_next)
2938  {
2939  ASSERT (c != super_c) ;
2940  ASSERT (COL_IS_ALIVE (c)) ;
2941  ASSERT (Col [c].shared3.hash == hash) ;
2942 
2943  /* not identical if lengths or scores are different */
2944  if (Col [c].length != length ||
2945  Col [c].shared2.score != Col [super_c].shared2.score)
2946  {
2947  prev_c = c ;
2948  continue ;
2949  }
2950 
2951  /* compare the two columns */
2952  cp1 = &A [Col [super_c].start] ;
2953  cp2 = &A [Col [c].start] ;
2954 
2955  for (i = 0 ; i < length ; i++)
2956  {
2957  /* the columns are "clean" (no dead rows) */
2958  ASSERT (ROW_IS_ALIVE (*cp1)) ;
2959  ASSERT (ROW_IS_ALIVE (*cp2)) ;
2960  /* row indices will same order for both supercols, */
2961  /* no gather scatter nessasary */
2962  if (*cp1++ != *cp2++)
2963  {
2964  break ;
2965  }
2966  }
2967 
2968  /* the two columns are different if the for-loop "broke" */
2969  if (i != length)
2970  {
2971  prev_c = c ;
2972  continue ;
2973  }
2974 
2975  /* === Got it! two columns are identical =================== */
2976 
2977  ASSERT (Col [c].shared2.score == Col [super_c].shared2.score) ;
2978 
2979  Col [super_c].shared1.thickness += Col [c].shared1.thickness ;
2980  Col [c].shared1.parent = super_c ;
2982  /* order c later, in order_children() */
2983  Col [c].shared2.order = EMPTY ;
2984  /* remove c from hash bucket */
2985  Col [prev_c].shared4.hash_next = Col [c].shared4.hash_next ;
2986  }
2987  }
2988 
2989  /* === Empty this hash bucket ======================================= */
2990 
2991  if (head_column > EMPTY)
2992  {
2993  /* corresponding degree list "hash" is not empty */
2994  Col [head_column].shared3.headhash = EMPTY ;
2995  }
2996  else
2997  {
2998  /* corresponding degree list "hash" is empty */
2999  head [hash] = EMPTY ;
3000  }
3001  }
3002 }
3003 
3004 
3005 /* ========================================================================== */
3006 /* === garbage_collection =================================================== */
3007 /* ========================================================================== */
3008 
3009 /*
3010  Defragments and compacts columns and rows in the workspace A. Used when
3011  all avaliable memory has been used while performing row merging. Returns
3012  the index of the first free position in A, after garbage collection. The
3013  time taken by this routine is linear is the size of the array A, which is
3014  itself linear in the number of nonzeros in the input matrix.
3015  Not user-callable.
3016 */
3017 
3018 PRIVATE Int garbage_collection /* returns the new value of pfree */
3020  /* === Parameters ======================================================= */
3021 
3022  Int n_row, /* number of rows */
3023  Int n_col, /* number of columns */
3024  Colamd_Row Row [], /* row info */
3025  Colamd_Col Col [], /* column info */
3026  Int A [], /* A [0 ... Alen-1] holds the matrix */
3027  Int *pfree /* &A [0] ... pfree is in use */
3028 )
3029 {
3030  /* === Local variables ================================================== */
3031 
3032  Int *psrc ; /* source pointer */
3033  Int *pdest ; /* destination pointer */
3034  Int j ; /* counter */
3035  Int r ; /* a row index */
3036  Int c ; /* a column index */
3037  Int length ; /* length of a row or column */
3038 
3039 #ifndef NDEBUG
3040  Int debug_rows ;
3041  DEBUG2 (("Defrag..\n")) ;
3042  for (psrc = &A[0] ; psrc < pfree ; psrc++) ASSERT (*psrc >= 0) ;
3043  debug_rows = 0 ;
3044 #endif /* NDEBUG */
3045 
3046  /* === Defragment the columns =========================================== */
3047 
3048  pdest = &A[0] ;
3049  for (c = 0 ; c < n_col ; c++)
3050  {
3051  if (COL_IS_ALIVE (c))
3052  {
3053  psrc = &A [Col [c].start] ;
3054 
3055  /* move and compact the column */
3056  ASSERT (pdest <= psrc) ;
3057  Col [c].start = (Int) (pdest - &A [0]) ;
3058  length = Col [c].length ;
3059  for (j = 0 ; j < length ; j++)
3060  {
3061  r = *psrc++ ;
3062  if (ROW_IS_ALIVE (r))
3063  {
3064  *pdest++ = r ;
3065  }
3066  }
3067  Col [c].length = (Int) (pdest - &A [Col [c].start]) ;
3068  }
3069  }
3070 
3071  /* === Prepare to defragment the rows =================================== */
3072 
3073  for (r = 0 ; r < n_row ; r++)
3074  {
3075  if (ROW_IS_DEAD (r) || (Row [r].length == 0))
3076  {
3077  /* This row is already dead, or is of zero length. Cannot compact
3078  * a row of zero length, so kill it. NOTE: in the current version,
3079  * there are no zero-length live rows. Kill the row (for the first
3080  * time, or again) just to be safe. */
3081  KILL_ROW (r) ;
3082  }
3083  else
3084  {
3085  /* save first column index in Row [r].shared2.first_column */
3086  psrc = &A [Row [r].start] ;
3087  Row [r].shared2.first_column = *psrc ;
3088  ASSERT (ROW_IS_ALIVE (r)) ;
3089  /* flag the start of the row with the one's complement of row */
3090  *psrc = ONES_COMPLEMENT (r) ;
3091 #ifndef NDEBUG
3092  debug_rows++ ;
3093 #endif /* NDEBUG */
3094  }
3095  }
3096 
3097  /* === Defragment the rows ============================================== */
3098 
3099  psrc = pdest ;
3100  while (psrc < pfree)
3101  {
3102  /* find a negative number ... the start of a row */
3103  if (*psrc++ < 0)
3104  {
3105  psrc-- ;
3106  /* get the row index */
3107  r = ONES_COMPLEMENT (*psrc) ;
3108  ASSERT (r >= 0 && r < n_row) ;
3109  /* restore first column index */
3110  *psrc = Row [r].shared2.first_column ;
3111  ASSERT (ROW_IS_ALIVE (r)) ;
3112  ASSERT (Row [r].length > 0) ;
3113  /* move and compact the row */
3114  ASSERT (pdest <= psrc) ;
3115  Row [r].start = (Int) (pdest - &A [0]) ;
3116  length = Row [r].length ;
3117  for (j = 0 ; j < length ; j++)
3118  {
3119  c = *psrc++ ;
3120  if (COL_IS_ALIVE (c))
3121  {
3122  *pdest++ = c ;
3123  }
3124  }
3125  Row [r].length = (Int) (pdest - &A [Row [r].start]) ;
3126  ASSERT (Row [r].length > 0) ;
3127 #ifndef NDEBUG
3128  debug_rows-- ;
3129 #endif /* NDEBUG */
3130  }
3131  }
3132  /* ensure we found all the rows */
3133  ASSERT (debug_rows == 0) ;
3134 
3135  /* === Return the new value of pfree ==================================== */
3136 
3137  return ((Int) (pdest - &A [0])) ;
3138 }
3139 
3140 
3141 /* ========================================================================== */
3142 /* === clear_mark =========================================================== */
3143 /* ========================================================================== */
3144 
3145 /*
3146  Clears the Row [].shared2.mark array, and returns the new tag_mark.
3147  Return value is the new tag_mark. Not user-callable.
3148 */
3149 
3150 PRIVATE Int clear_mark /* return the new value for tag_mark */
3152  /* === Parameters ======================================================= */
3153 
3154  Int tag_mark, /* new value of tag_mark */
3155  Int max_mark, /* max allowed value of tag_mark */
3156 
3157  Int n_row, /* number of rows in A */
3158  Colamd_Row Row [] /* Row [0 ... n_row-1].shared2.mark is set to zero */
3159 )
3160 {
3161  /* === Local variables ================================================== */
3162 
3163  Int r ;
3164 
3165  if (tag_mark <= 0 || tag_mark >= max_mark)
3166  {
3167  for (r = 0 ; r < n_row ; r++)
3168  {
3169  if (ROW_IS_ALIVE (r))
3170  {
3171  Row [r].shared2.mark = 0 ;
3172  }
3173  }
3174  tag_mark = 1 ;
3175  }
3176 
3177  return (tag_mark) ;
3178 }
3179 
3180 
3181 /* ========================================================================== */
3182 /* === print_report ========================================================= */
3183 /* ========================================================================== */
3184 
3185 PRIVATE void print_report
3187  char *method,
3188  Int stats [COLAMD_STATS]
3189 )
3190 {
3191 
3192  Int i1, i2, i3 ;
3193 
3194  PRINTF (("\n%s version %d.%d, %s: ", method,
3196 
3197  if (!stats)
3198  {
3199  PRINTF (("No statistics available.\n")) ;
3200  return ;
3201  }
3202 
3203  i1 = stats [COLAMD_INFO1] ;
3204  i2 = stats [COLAMD_INFO2] ;
3205  i3 = stats [COLAMD_INFO3] ;
3206 
3207  if (stats [COLAMD_STATUS] >= 0)
3208  {
3209  PRINTF (("OK. ")) ;
3210  }
3211  else
3212  {
3213  PRINTF (("ERROR. ")) ;
3214  }
3215 
3216  switch (stats [COLAMD_STATUS])
3217  {
3218 
3219  case COLAMD_OK_BUT_JUMBLED:
3220 
3221  PRINTF(("Matrix has unsorted or duplicate row indices.\n")) ;
3222 
3223  PRINTF(("%s: number of duplicate or out-of-order row indices: %d\n",
3224  method, i3)) ;
3225 
3226  PRINTF(("%s: last seen duplicate or out-of-order row index: %d\n",
3227  method, INDEX (i2))) ;
3228 
3229  PRINTF(("%s: last seen in column: %d",
3230  method, INDEX (i1))) ;
3231 
3232  /* no break - fall through to next case instead */
3233 
3234  case COLAMD_OK:
3235 
3236  PRINTF(("\n")) ;
3237 
3238  PRINTF(("%s: number of dense or empty rows ignored: %d\n",
3239  method, stats [COLAMD_DENSE_ROW])) ;
3240 
3241  PRINTF(("%s: number of dense or empty columns ignored: %d\n",
3242  method, stats [COLAMD_DENSE_COL])) ;
3243 
3244  PRINTF(("%s: number of garbage collections performed: %d\n",
3245  method, stats [COLAMD_DEFRAG_COUNT])) ;
3246  break ;
3247 
3249 
3250  PRINTF(("Array A (row indices of matrix) not present.\n")) ;
3251  break ;
3252 
3254 
3255  PRINTF(("Array p (column pointers for matrix) not present.\n")) ;
3256  break ;
3257 
3259 
3260  PRINTF(("Invalid number of rows (%d).\n", i1)) ;
3261  break ;
3262 
3264 
3265  PRINTF(("Invalid number of columns (%d).\n", i1)) ;
3266  break ;
3267 
3269 
3270  PRINTF(("Invalid number of nonzero entries (%d).\n", i1)) ;
3271  break ;
3272 
3274 
3275  PRINTF(("Invalid column pointer, p [0] = %d, must be zero.\n", i1));
3276  break ;
3277 
3279 
3280  PRINTF(("Array A too small.\n")) ;
3281  PRINTF((" Need Alen >= %d, but given only Alen = %d.\n",
3282  i1, i2)) ;
3283  break ;
3284 
3286 
3287  PRINTF
3288  (("Column %d has a negative number of nonzero entries (%d).\n",
3289  INDEX (i1), i2)) ;
3290  break ;
3291 
3293 
3294  PRINTF
3295  (("Row index (row %d) out of bounds (%d to %d) in column %d.\n",
3296  INDEX (i2), INDEX (0), INDEX (i3-1), INDEX (i1))) ;
3297  break ;
3298 
3300 
3301  PRINTF(("Out of memory.\n")) ;
3302  break ;
3303 
3304  /* v2.4: internal-error case deleted */
3305  }
3306 }
3307 
3308 
3309 
3310 
3311 /* ========================================================================== */
3312 /* === colamd debugging routines ============================================ */
3313 /* ========================================================================== */
3314 
3315 /* When debugging is disabled, the remainder of this file is ignored. */
3316 
3317 #ifndef NDEBUG
3318 
3319 
3320 /* ========================================================================== */
3321 /* === debug_structures ===================================================== */
3322 /* ========================================================================== */
3323 
3324 /*
3325  At this point, all empty rows and columns are dead. All live columns
3326  are "clean" (containing no dead rows) and simplicial (no supercolumns
3327  yet). Rows may contain dead columns, but all live rows contain at
3328  least one live column.
3329 */
3330 
3331 PRIVATE void debug_structures
3332 (
3333  /* === Parameters ======================================================= */
3334 
3335  Int n_row,
3336  Int n_col,
3337  Colamd_Row Row [],
3338  Colamd_Col Col [],
3339  Int A [],
3340  Int n_col2
3341 )
3342 {
3343  /* === Local variables ================================================== */
3344 
3345  Int i ;
3346  Int c ;
3347  Int *cp ;
3348  Int *cp_end ;
3349  Int len ;
3350  Int score ;
3351  Int r ;
3352  Int *rp ;
3353  Int *rp_end ;
3354  Int deg ;
3355 
3356  /* === Check A, Row, and Col ============================================ */
3357 
3358  for (c = 0 ; c < n_col ; c++)
3359  {
3360  if (COL_IS_ALIVE (c))
3361  {
3362  len = Col [c].length ;
3363  score = Col [c].shared2.score ;
3364  DEBUG4 (("initial live col %5d %5d %5d\n", c, len, score)) ;
3365  ASSERT (len > 0) ;
3366  ASSERT (score >= 0) ;
3367  ASSERT (Col [c].shared1.thickness == 1) ;
3368  cp = &A [Col [c].start] ;
3369  cp_end = cp + len ;
3370  while (cp < cp_end)
3371  {
3372  r = *cp++ ;
3373  ASSERT (ROW_IS_ALIVE (r)) ;
3374  }
3375  }
3376  else
3377  {
3378  i = Col [c].shared2.order ;
3379  ASSERT (i >= n_col2 && i < n_col) ;
3380  }
3381  }
3382 
3383  for (r = 0 ; r < n_row ; r++)
3384  {
3385  if (ROW_IS_ALIVE (r))
3386  {
3387  i = 0 ;
3388  len = Row [r].length ;
3389  deg = Row [r].shared1.degree ;
3390  ASSERT (len > 0) ;
3391  ASSERT (deg > 0) ;
3392  rp = &A [Row [r].start] ;
3393  rp_end = rp + len ;
3394  while (rp < rp_end)
3395  {
3396  c = *rp++ ;
3397  if (COL_IS_ALIVE (c))
3398  {
3399  i++ ;
3400  }
3401  }
3402  ASSERT (i > 0) ;
3403  }
3404  }
3405 }
3406 
3407 
3408 /* ========================================================================== */
3409 /* === debug_deg_lists ====================================================== */
3410 /* ========================================================================== */
3411 
3412 /*
3413  Prints the contents of the degree lists. Counts the number of columns
3414  in the degree list and compares it to the total it should have. Also
3415  checks the row degrees.
3416 */
3417 
3418 PRIVATE void debug_deg_lists
3419 (
3420  /* === Parameters ======================================================= */
3421 
3422  Int n_row,
3423  Int n_col,
3424  Colamd_Row Row [],
3425  Colamd_Col Col [],
3426  Int head [],
3427  Int min_score,
3428  Int should,
3429  Int max_deg
3430 )
3431 {
3432  /* === Local variables ================================================== */
3433 
3434  Int deg ;
3435  Int col ;
3436  Int have ;
3437  Int row ;
3438 
3439  /* === Check the degree lists =========================================== */
3440 
3441  if (n_col > 10000 && colamd_debug <= 0)
3442  {
3443  return ;
3444  }
3445  have = 0 ;
3446  DEBUG4 (("Degree lists: %d\n", min_score)) ;
3447  for (deg = 0 ; deg <= n_col ; deg++)
3448  {
3449  col = head [deg] ;
3450  if (col == EMPTY)
3451  {
3452  continue ;
3453  }
3454  DEBUG4 (("%d:", deg)) ;
3455  while (col != EMPTY)
3456  {
3457  DEBUG4 ((" %d", col)) ;
3458  have += Col [col].shared1.thickness ;
3459  ASSERT (COL_IS_ALIVE (col)) ;
3460  col = Col [col].shared4.degree_next ;
3461  }
3462  DEBUG4 (("\n")) ;
3463  }
3464  DEBUG4 (("should %d have %d\n", should, have)) ;
3465  ASSERT (should == have) ;
3466 
3467  /* === Check the row degrees ============================================ */
3468 
3469  if (n_row > 10000 && colamd_debug <= 0)
3470  {
3471  return ;
3472  }
3473  for (row = 0 ; row < n_row ; row++)
3474  {
3475  if (ROW_IS_ALIVE (row))
3476  {
3477  ASSERT (Row [row].shared1.degree <= max_deg) ;
3478  }
3479  }
3480 }
3481 
3482 
3483 /* ========================================================================== */
3484 /* === debug_mark =========================================================== */
3485 /* ========================================================================== */
3486 
3487 /*
3488  Ensures that the tag_mark is less that the maximum and also ensures that
3489  each entry in the mark array is less than the tag mark.
3490 */
3491 
3492 PRIVATE void debug_mark
3493 (
3494  /* === Parameters ======================================================= */
3495 
3496  Int n_row,
3497  Colamd_Row Row [],
3498  Int tag_mark,
3499  Int max_mark
3500 )
3501 {
3502  /* === Local variables ================================================== */
3503 
3504  Int r ;
3505 
3506  /* === Check the Row marks ============================================== */
3507 
3508  ASSERT (tag_mark > 0 && tag_mark <= max_mark) ;
3509  if (n_row > 10000 && colamd_debug <= 0)
3510  {
3511  return ;
3512  }
3513  for (r = 0 ; r < n_row ; r++)
3514  {
3515  ASSERT (Row [r].shared2.mark < tag_mark) ;
3516  }
3517 }
3518 
3519 
3520 /* ========================================================================== */
3521 /* === debug_matrix ========================================================= */
3522 /* ========================================================================== */
3523 
3524 /*
3525  Prints out the contents of the columns and the rows.
3526 */
3527 
3528 PRIVATE void debug_matrix
3529 (
3530  /* === Parameters ======================================================= */
3531 
3532  Int n_row,
3533  Int n_col,
3534  Colamd_Row Row [],
3535  Colamd_Col Col [],
3536  Int A []
3537 )
3538 {
3539  /* === Local variables ================================================== */
3540 
3541  Int r ;
3542  Int c ;
3543  Int *rp ;
3544  Int *rp_end ;
3545  Int *cp ;
3546  Int *cp_end ;
3547 
3548  /* === Dump the rows and columns of the matrix ========================== */
3549 
3550  if (colamd_debug < 3)
3551  {
3552  return ;
3553  }
3554  DEBUG3 (("DUMP MATRIX:\n")) ;
3555  for (r = 0 ; r < n_row ; r++)
3556  {
3557  DEBUG3 (("Row %d alive? %d\n", r, ROW_IS_ALIVE (r))) ;
3558  if (ROW_IS_DEAD (r))
3559  {
3560  continue ;
3561  }
3562  DEBUG3 (("start %d length %d degree %d\n",
3563  Row [r].start, Row [r].length, Row [r].shared1.degree)) ;
3564  rp = &A [Row [r].start] ;
3565  rp_end = rp + Row [r].length ;
3566  while (rp < rp_end)
3567  {
3568  c = *rp++ ;
3569  DEBUG4 ((" %d col %d\n", COL_IS_ALIVE (c), c)) ;
3570  }
3571  }
3572 
3573  for (c = 0 ; c < n_col ; c++)
3574  {
3575  DEBUG3 (("Col %d alive? %d\n", c, COL_IS_ALIVE (c))) ;
3576  if (COL_IS_DEAD (c))
3577  {
3578  continue ;
3579  }
3580  DEBUG3 (("start %d length %d shared1 %d shared2 %d\n",
3581  Col [c].start, Col [c].length,
3582  Col [c].shared1.thickness, Col [c].shared2.score)) ;
3583  cp = &A [Col [c].start] ;
3584  cp_end = cp + Col [c].length ;
3585  while (cp < cp_end)
3586  {
3587  r = *cp++ ;
3588  DEBUG4 ((" %d row %d\n", ROW_IS_ALIVE (r), r)) ;
3589  }
3590  }
3591 }
3592 
3593 PRIVATE void colamd_get_debug
3594 (
3595  char *method
3596 )
3597 {
3598  FILE *f ;
3599  colamd_debug = 0 ; /* no debug printing */
3600  f = fopen ("debug", "r") ;
3601  if (f == (FILE *) NULL)
3602  {
3603  colamd_debug = 0 ;
3604  }
3605  else
3606  {
3607  fscanf (f, "%d", &colamd_debug) ;
3608  fclose (f) ;
3609  }
3610  DEBUG0 (("%s: debug version, D = %d (THIS WILL BE SLOW!)\n",
3611  method, colamd_debug)) ;
3612 }
3613 
3614 #endif /* NDEBUG */
#define DEBUG3(params)
union Colamd_Row_struct::@17 shared2
PRIVATE void print_report(char *method, Int stats[COLAMD_STATS])
#define COLAMD_INFO1
#define DEBUG2(params)
void f()
#define COLAMD_ERROR_ncol_negative
#define COLAMD_MAIN_VERSION
Definition: amesos_colamd.h:86
UF_long CHOLMOD() nnz(cholmod_sparse *A, cholmod_common *Common)
static size_t t_add(size_t a, size_t b, int *ok)
#define COLAMD_ERROR_p0_nonzero
#define COLAMD_KNOBS
Definition: amesos_colamd.h:97
#define COL_IS_DEAD(c)
#define COL_IS_DEAD_PRINCIPAL(c)
#define ROW_IS_DEAD(r)
#define COLAMD_ERROR_row_index_out_of_bounds
#define Int_MAX
#define COLAMD_INFO2
union Colamd_Row_struct::@16 shared1
PRIVATE Int init_rows_cols(Int n_row, Int n_col, Colamd_Row Row[], Colamd_Col Col[], Int A[], Int p[], Int stats[COLAMD_STATS])
PRIVATE Int garbage_collection(Int n_row, Int n_col, Colamd_Row Row[], Colamd_Col Col[], Int A[], Int *pfree)
#define COLAMD_STATUS
struct Colamd_Col_struct Colamd_Col
#define MIN(a, b)
#define ROW_IS_ALIVE(r)
#define COLAMD_ERROR_nrow_negative
#define COLAMD_ERROR_out_of_memory
#define SYMAMD_MAIN
#define COLAMD_OK_BUT_JUMBLED
#define ROW_IS_MARKED_DEAD(row_mark)
#define ASSERT(expression)
#define COLAMD_INFO3
#define COLAMD_C(n_col, ok)
PRIVATE void order_children(Int n_col, Colamd_Col Col[], Int p[])
#define COLAMD_set_defaults
#define COLAMD_AGGRESSIVE
#define COLAMD_R(n_row, ok)
#define COLAMD_DEFRAG_COUNT
#define NDEBUG
#define COLAMD_MAIN
#define KILL_ROW(r)
#define INDEX(i)
union Colamd_Col_struct::@13 shared2
#define KILL_PRINCIPAL_COL(c)
#define COLAMD_ERROR_nnz_negative
PRIVATE void detect_super_cols(Colamd_Col Col[], Int A[], Int head[], Int row_start, Int row_length)
int CHOLMOD() start(cholmod_common *Common)
#define COLAMD_ERROR_A_not_present
#define COLAMD_ERROR_A_too_small
#define COLAMD_recommended
#define COLAMD_report
#define FALSE
#define NULL
#define COLAMD_OK
union Colamd_Col_struct::@12 shared1
#define COLAMD_STATS
PRIVATE Int find_ordering(Int n_row, Int n_col, Int Alen, Colamd_Row Row[], Colamd_Col Col[], Int A[], Int head[], Int n_col2, Int max_deg, Int pfree, Int aggressive)
#define COLAMD_SUB_VERSION
Definition: amesos_colamd.h:87
#define PRIVATE
struct Colamd_Row_struct Colamd_Row
#define TRUE
#define DEBUG4(params)
#define DEBUG1(params)
#define PRINTF(params)
PRIVATE Int clear_mark(Int tag_mark, Int max_mark, Int n_row, Colamd_Row Row[])
#define COLAMD_DENSE_ROW
#define COLAMD_DATE
Definition: amesos_colamd.h:84
#define DEBUG0(params)
#define COLAMD_ERROR_col_length_negative
#define KILL_NON_PRINCIPAL_COL(c)
#define EMPTY
#define COL_IS_ALIVE(c)
static size_t t_mult(size_t a, size_t k, int *ok)
#define MAX(a, b)
#define DENSE_DEGREE(alpha, n)
PRIVATE void init_scoring(Int n_row, Int n_col, Colamd_Row Row[], Colamd_Col Col[], Int A[], Int head[], double knobs[COLAMD_KNOBS], Int *p_n_row2, Int *p_n_col2, Int *p_max_deg)
union Colamd_Col_struct::@15 shared4
union Colamd_Col_struct::@14 shared3
int n
#define SYMAMD_report
#define Int
#define COLAMD_DENSE_COL
#define PUBLIC
#define COLAMD_ERROR_p_not_present
#define ONES_COMPLEMENT(r)