Fixing pointers and debug formatting.
[COLAMDJ.git] / src / main / java / edu / ufl / cise / colamd / tdouble / Dcolamd.java
blobc6035890cf3e06fa7cb92f23553893b614d18bdd
1 /**
2 * Copyright (c) 1998-2007, Timothy A. Davis, All Rights Reserved.
4 * This library is free software; you can redistribute it and/or
5 * modify it under the terms of the GNU Lesser General Public
6 * License as published by the Free Software Foundation; either
7 * version 2.1 of the License, or (at your option) any later version.
9 * This library is distributed in the hope that it will be useful,
10 * but WITHOUT ANY WARRANTY; without even the implied warranty of
11 * MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE. See the GNU
12 * Lesser General Public License for more details.
14 * You should have received a copy of the GNU Lesser General Public
15 * License along with this library; if not, write to the Free Software
16 * Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA 02110-1301
19 package edu.ufl.cise.colamd.tdouble;
21 import java.io.BufferedReader;
22 import java.io.File;
23 import java.io.FileReader;
24 import java.io.IOException;
26 /**
27 * COLAMD / SYMAMD
29 * colamd: an approximate minimum degree column ordering algorithm,
30 * for LU factorization of symmetric or unsymmetric matrices,
31 * QR factorization, least squares, interior point methods for
32 * linear programming problems, and other related problems.
34 * symamd: an approximate minimum degree ordering algorithm for Cholesky
35 * factorization of symmetric matrices.
37 * Purpose:
39 * Colamd computes a permutation Q such that the Cholesky factorization of
40 * (AQ)'(AQ) has less fill-in and requires fewer floating point operations
41 * than A'A. This also provides a good ordering for sparse partial
42 * pivoting methods, P(AQ) = LU, where Q is computed prior to numerical
43 * factorization, and P is computed during numerical factorization via
44 * conventional partial pivoting with row interchanges. Colamd is the
45 * column ordering method used in SuperLU, part of the ScaLAPACK library.
46 * It is also available as built-in function in MATLAB Version 6,
47 * available from MathWorks, Inc. (http://www.mathworks.com). This
48 * routine can be used in place of colmmd in MATLAB.
50 * Symamd computes a permutation P of a symmetric matrix A such that the
51 * Cholesky factorization of PAP' has less fill-in and requires fewer
52 * floating point operations than A. Symamd constructs a matrix M such
53 * that M'M has the same nonzero pattern of A, and then orders the columns
54 * of M using colmmd. The column ordering of M is then returned as the
55 * row and column ordering P of A.
57 * Authors:
59 * The authors of the code itself are Stefan I. Larimore and Timothy A.
60 * Davis (davis at cise.ufl.edu), University of Florida. The algorithm was
61 * developed in collaboration with John Gilbert, Xerox PARC, and Esmond
62 * Ng, Oak Ridge National Laboratory.
64 * Acknowledgements:
66 * This work was supported by the National Science Foundation, under
67 * grants DMS-9504974 and DMS-9803599.
69 public class Dcolamd {
71 /* ========================================================================== */
72 /* === COLAMD version ======================================================= */
73 /* ========================================================================== */
75 /* COLAMD Version 2.4 and later will include the following definitions.
76 * As an example, to test if the version you are using is 2.4 or later:
78 * if (COLAMD_VERSION >= COLAMD_VERSION_CODE (2,4)) ...
81 public static String COLAMD_DATE = "Jan 25, 2011" ;
82 public static int COLAMD_VERSION_CODE(int main, int sub)
84 return ((main) * 1000 + (sub)) ;
86 public static int COLAMD_MAIN_VERSION = 2 ;
87 public static int COLAMD_SUB_VERSION = 7 ;
88 public static int COLAMD_SUBSUB_VERSION = 3 ;
89 public static int COLAMD_VERSION = COLAMD_VERSION_CODE(COLAMD_MAIN_VERSION,
90 COLAMD_SUB_VERSION) ;
92 /* ========================================================================== */
93 /* === Knob and statistics definitions ====================================== */
94 /* ========================================================================== */
96 /* size of the knobs [ ] array. Only knobs [0..1] are currently used. */
97 public static int COLAMD_KNOBS = 20 ;
99 /* number of output statistics. Only stats [0..6] are currently used. */
100 public static int COLAMD_STATS = 20 ;
102 /* knobs [0] and stats [0]: dense row knob and output statistic. */
103 public static int COLAMD_DENSE_ROW = 0 ;
105 /* knobs [1] and stats [1]: dense column knob and output statistic. */
106 public static int COLAMD_DENSE_COL = 1 ;
108 /* knobs [2]: aggressive absorption */
109 public static int COLAMD_AGGRESSIVE = 2 ;
111 /* stats [2]: memory defragmentation count output statistic */
112 public static int COLAMD_DEFRAG_COUNT = 2 ;
114 /* stats [3]: colamd status: zero OK, > 0 warning or notice, < 0 error */
115 public static int COLAMD_STATUS = 3 ;
117 /* stats [4..6]: error info, or info on jumbled columns */
118 public static final int COLAMD_INFO1 = 4 ;
119 public static final int COLAMD_INFO2 = 5 ;
120 public static final int COLAMD_INFO3 = 6 ;
122 /* error codes returned in stats [3]: */
123 public static final int COLAMD_OK = (0) ;
124 public static final int COLAMD_OK_BUT_JUMBLED = (1) ;
125 public static final int COLAMD_ERROR_A_not_present = (-1) ;
126 public static final int COLAMD_ERROR_p_not_present = (-2) ;
127 public static final int COLAMD_ERROR_nrow_negative = (-3) ;
128 public static final int COLAMD_ERROR_ncol_negative = (-4) ;
129 public static final int COLAMD_ERROR_nnz_negative = (-5) ;
130 public static final int COLAMD_ERROR_p0_nonzero = (-6) ;
131 public static final int COLAMD_ERROR_A_too_small = (-7) ;
132 public static final int COLAMD_ERROR_col_length_negative = (-8) ;
133 public static final int COLAMD_ERROR_row_index_out_of_bounds = (-9) ;
134 public static final int COLAMD_ERROR_out_of_memory = (-10) ;
135 public static final int COLAMD_ERROR_internal_error = (-999) ;
137 /* ========================================================================== */
138 /* === Scaffolding code definitions ======================================== */
139 /* ========================================================================== */
141 public static boolean NDEBUG = true;
143 public static boolean NPRINT = true;
146 Our "scaffolding code" philosophy: In our opinion, well-written library
147 code should keep its "debugging" code, and just normally have it turned off
148 by the compiler so as not to interfere with performance. This serves
149 several purposes:
151 (1) assertions act as comments to the reader, telling you what the code
152 expects at that point. All assertions will always be true (unless
153 there really is a bug, of course).
155 (2) leaving in the scaffolding code assists anyone who would like to modify
156 the code, or understand the algorithm (by reading the debugging output,
157 one can get a glimpse into what the code is doing).
159 (3) (gasp!) for actually finding bugs. This code has been heavily tested
160 and "should" be fully functional and bug-free ... but you never know...
162 The code will become outrageously slow when debugging is
163 enabled. To control the level of debugging output, set an environment
164 variable D to 0 (little), 1 (some), 2, 3, or 4 (lots). When debugging,
165 you should see the following message on the standard output:
167 colamd: debug version, D = 1 (THIS WILL BE SLOW!)
169 or a similar message for symamd. If you don't, then debugging has not
170 been enabled.
174 private static int Int_MAX = Integer.MAX_VALUE;
176 /* ========================================================================== */
177 /* === Definitions ========================================================== */
178 /* ========================================================================== */
180 protected static double sqrt (double a)
182 return Math.sqrt (a) ;
185 private static final int MAX (int a, int b)
187 return (((a) > (b)) ? (a) : (b)) ;
191 private static final int MIN (int a, int b)
193 return (((a) < (b)) ? (a) : (b)) ;
196 private static final double MAX (double a, double b)
198 return (((a) > (b)) ? (a) : (b)) ;
202 private static int DENSE_DEGREE (double alpha, int n)
204 return ((int) MAX (16.0, (alpha) * sqrt ((double) (n)))) ;
207 private static int ONES_COMPLEMENT (int r)
209 return (-(r)-1) ;
212 private static final int TRUE = (1) ;
213 private static final int FALSE = (0) ;
215 private static final int EMPTY = (-1) ;
217 /* Row and column status */
218 private static final int ALIVE = (0) ;
219 private static final int DEAD = (-1) ;
221 /* Column status */
222 private static final int DEAD_PRINCIPAL = (-1) ;
223 private static final int DEAD_NON_PRINCIPAL = (-2) ;
225 /* Row and column status update and checking. */
226 private static boolean ROW_IS_DEAD(Colamd_Row[] Row, int r)
228 return ROW_IS_MARKED_DEAD (Row [r].mark) ;
230 private static boolean ROW_IS_MARKED_DEAD(int row_mark)
232 return (row_mark < ALIVE) ;
234 private static boolean ROW_IS_ALIVE(Colamd_Row[] Row, int r)
236 return (Row [r].mark >= ALIVE) ;
238 private static boolean COL_IS_DEAD(Colamd_Col[] Col, int c)
240 return (Col [c].start < ALIVE) ;
242 private static boolean COL_IS_ALIVE(Colamd_Col[] Col, int c)
244 return (Col [c].start >= ALIVE) ;
246 private static boolean COL_IS_DEAD_PRINCIPAL(Colamd_Col[] Col, int c)
248 return (Col [c].start == DEAD_PRINCIPAL) ;
250 private static void KILL_ROW(Colamd_Row[] Row, int r)
252 Row [r].mark = DEAD ;
254 private static void KILL_PRINCIPAL_COL(Colamd_Col[] Col, int c)
256 Col [c].start = DEAD_PRINCIPAL ;
258 private static void KILL_NON_PRINCIPAL_COL(Colamd_Col[] Col, int c)
260 Col [c].start = DEAD_NON_PRINCIPAL ;
263 /* ========================================================================== */
264 /* === Colamd reporting mechanism =========================================== */
265 /* ========================================================================== */
268 * In Java, matrices are 0-based and indices are reported as such
269 * in *_report.
271 private static int INDEX (int i)
273 return (i) ;
277 * All output goes through PRINTF.
279 private static void PRINTF (String format, Object... args)
281 if (!NPRINT)
283 System.out.printf (format, args) ;
287 /** debug print level */
288 public static int colamd_debug = 0 ;
290 protected static void DEBUG0 (String format, Object... args)
292 if (!NDEBUG)
294 PRINTF (format, args) ;
298 protected static void DEBUG1(String format, Object... args)
300 if (!NDEBUG)
302 if (colamd_debug >= 1) PRINTF (format, args) ;
306 protected static void DEBUG2(String format, Object... args)
308 if (!NDEBUG)
310 if (colamd_debug >= 2) PRINTF (format, args) ;
314 protected static void DEBUG3(String format, Object... args)
316 if (!NDEBUG)
318 if (colamd_debug >= 3) PRINTF (format, args) ;
322 protected static void DEBUG4(String format, Object... args)
324 if (!NDEBUG)
326 if (colamd_debug >= 4) PRINTF (format, args) ;
330 protected static void ASSERT (boolean a)
332 if (!NDEBUG)
334 assert a ;
338 protected static void ASSERT (int a)
340 ASSERT (a != 0) ;
343 /* ========================================================================== */
344 /* === USER-CALLABLE ROUTINES: ============================================== */
345 /* ========================================================================== */
347 /* ========================================================================== */
348 /* === colamd_recommended =================================================== */
349 /* ========================================================================== */
352 The colamd_recommended routine returns the suggested size for Alen. This
353 value has been determined to provide good balance between the number of
354 garbage collections and the memory requirements for colamd. If any
355 argument is negative, or if integer overflow occurs, a 0 is returned as an
356 error condition. 2*nnz space is required for the row and column
357 indices of the matrix. COLAMD_C (n_col) + COLAMD_R (n_row) space is
358 required for the Col and Row arrays, respectively, which are internal to
359 colamd (roughly 6*n_col + 4*n_row). An additional n_col space is the
360 minimal amount of "elbow room", and nnz/5 more space is recommended for
361 run time efficiency.
363 Alen is approximately 2.2*nnz + 7*n_col + 4*n_row + 10.
365 This function is not needed when using symamd.
369 * add two values of type int, and check for integer overflow
371 private static int t_add (int a, int b, int[] ok)
373 (ok[0]) = (ok[0] != 0) && ((a + b) >= MAX (a,b)) ? 1 : 0;
374 return ((ok[0] != 0) ? (a + b) : 0) ;
378 * compute a*k where k is a small integer, and check for integer overflow
380 static int t_mult (int a, int k, int[] ok)
382 int i, s = 0 ;
383 for (i = 0 ; i < k ; i++)
385 s = t_add (s, a, ok) ;
387 return (s) ;
391 * size of the Col and Row structures
393 private static int COLAMD_C(int n_col, int[] ok)
395 // return ((t_mult (t_add (n_col, 1, ok), sizeof (Colamd_Col), ok) / sizeof (int))) ;
396 return t_add (n_col, 1, ok) ;
399 private static int COLAMD_R(int n_row, int[] ok)
401 // return ((t_mult (t_add (n_row, 1, ok), sizeof (Colamd_Row), ok) / sizeof (int))) ;
402 return t_add (n_row, 1, ok) ;
407 * @param nnz number of nonzeros in A. This must
408 * be the same value as p [n_col] in the call to
409 * colamd - otherwise you will get a wrong value
410 * of the recommended memory to use.
411 * @param n_row number of rows in A
412 * @param n_col number of columns in A
413 * @return recommended value of Alen. Returns 0
414 * if any input argument is negative. The use of this routine
415 * is optional. Not needed for symamd, which dynamically allocates
416 * its own memory.
418 public static int COLAMD_recommended (int nnz, int n_row, int n_col)
420 int s, c, r ;
421 int[] ok = new int [] { TRUE } ;
422 if (nnz < 0 || n_row < 0 || n_col < 0)
424 return (0) ;
426 s = t_mult (nnz, 2, ok) ; /* 2*nnz */
427 // c = COLAMD_C (n_col, ok) ; /* size of column structures */
428 // r = COLAMD_R (n_row, ok) ; /* size of row structures */
429 // s = t_add (s, c, ok) ;
430 // s = t_add (s, r, ok) ;
431 s = t_add (s, n_col, ok) ; /* elbow room */
432 s = t_add (s, nnz/5, ok) ; /* elbow room */
433 ok[0] = (s < Int_MAX) ? 1 : 0;
434 return (ok[0] != 0 ? s : 0) ;
438 /* ========================================================================== */
439 /* === colamd_set_defaults ================================================== */
440 /* ========================================================================== */
443 The colamd_set_defaults routine sets the default values of the user-
444 controllable parameters for colamd and symamd:
446 Colamd: rows with more than max (16, knobs [0] * sqrt (n_col))
447 entries are removed prior to ordering. Columns with more than
448 max (16, knobs [1] * sqrt (MIN (n_row,n_col))) entries are removed
449 prior to ordering, and placed last in the output column ordering.
451 Symamd: Rows and columns with more than max (16, knobs [0] * sqrt (n))
452 entries are removed prior to ordering, and placed last in the
453 output ordering.
455 knobs [0] dense row control
457 knobs [1] dense column control
459 knobs [2] if nonzero, do aggresive absorption
461 knobs [3..19] unused, but future versions might use this
466 * knobs [0] and knobs [1] control dense row and col detection:
468 * Colamd: rows with more than
469 * max (16, knobs [COLAMD_DENSE_ROW] * sqrt (n_col))
470 * entries are removed prior to ordering. Columns with more than
471 * max (16, knobs [COLAMD_DENSE_COL] * sqrt (MIN (n_row,n_col)))
472 * entries are removed prior to
473 * ordering, and placed last in the output column ordering.
475 * Symamd: uses only knobs [COLAMD_DENSE_ROW], which is knobs [0].
476 * Rows and columns with more than
477 * max (16, knobs [COLAMD_DENSE_ROW] * sqrt (n))
478 * entries are removed prior to ordering, and placed last in the
479 * output ordering.
481 * COLAMD_DENSE_ROW and COLAMD_DENSE_COL are defined as 0 and 1,
482 * respectively, in colamd.h. Default values of these two knobs
483 * are both 10. Currently, only knobs [0] and knobs [1] are
484 * used, but future versions may use more knobs. If so, they will
485 * be properly set to their defaults by the future version of
486 * colamd_set_defaults, so that the code that calls colamd will
487 * not need to change, assuming that you either use
488 * colamd_set_defaults, or pass a (double *) NULL pointer as the
489 * knobs array to colamd or symamd.
491 * knobs [2]: aggressive absorption
493 * knobs [COLAMD_AGGRESSIVE] controls whether or not to do
494 * aggressive absorption during the ordering. Default is TRUE.
496 * @param knobs knob array
498 public static void COLAMD_set_defaults (double[] knobs)
500 /* === Local variables ============================================== */
502 int i ;
504 if (knobs == null || knobs.length == 0)
506 return ; /* no knobs to initialize */
508 for (i = 0 ; i < COLAMD_KNOBS ; i++)
510 knobs [i] = 0 ;
512 knobs [COLAMD_DENSE_ROW] = 10 ;
513 knobs [COLAMD_DENSE_COL] = 10 ;
514 knobs [COLAMD_AGGRESSIVE] = TRUE ; /* default: do aggressive absorption*/
517 /* ========================================================================== */
518 /* === symamd =============================================================== */
519 /* ========================================================================== */
522 * The symamd routine computes an ordering P of a symmetric sparse
523 * matrix A such that the Cholesky factorization PAP' = LL' remains
524 * sparse. It is based on a column ordering of a matrix M constructed
525 * so that the nonzero pattern of M'M is the same as A. The matrix A
526 * is assumed to be symmetric; only the strictly lower triangular part
527 * is accessed. You must pass your selected memory allocator (usually
528 * calloc/free or mxCalloc/mxFree) to symamd, for it to allocate
529 * memory for the temporary matrix M.
531 * @param n number of rows and columns of A. Restriction: n >= 0.
532 * Symamd returns FALSE if n is negative.
533 * @param A an integer array of size nnz, where nnz = p [n].
535 * The row indices of the entries in column c of the matrix are
536 * held in A [(p [c]) ... (p [c+1]-1)]. The row indices in a
537 * given column c need not be in ascending order, and duplicate
538 * row indices may be present. However, symamd will run faster
539 * if the columns are in sorted order with no duplicate entries.
541 * The matrix is 0-based. That is, rows are in the range 0 to
542 * n-1, and columns are in the range 0 to n-1. Symamd
543 * returns FALSE if any row index is out of range.
545 * The contents of A are not modified.
546 * @param an integer array of size n+1. On input, it holds the
547 * "pointers" for the column form of the matrix A. Column c of
548 * the matrix A is held in A [(p [c]) ... (p [c+1]-1)]. The first
549 * entry, p [0], must be zero, and p [c] <= p [c+1] must hold
550 * for all c in the range 0 to n-1. The value p [n] is
551 * thus the total number of entries in the pattern of the matrix A.
552 * Symamd returns FALSE if these conditions are not met.
554 * The contents of p are not modified.
555 * @param perm On output, if symamd returns TRUE, the array perm holds the
556 * permutation P, where perm [0] is the first index in the new
557 * ordering, and perm [n-1] is the last. That is, perm [k] = j
558 * means that row and column j of A is the kth column in PAP',
559 * where k is in the range 0 to n-1 (perm [0] = j means
560 * that row and column j of A are the first row and column in
561 * PAP'). The array is used as a workspace during the ordering,
562 * which is why it must be of length n+1, not just n.
563 * @param knobs parameters (uses defaults if NULL)
564 * @param stats Statistics on the ordering, and error status.
565 * Symamd returns FALSE if stats is not present.
567 * stats [0]: number of dense or empty row and columns ignored
568 * (and ordered last in the output permutation
569 * perm). Note that a row/column can become
570 * "empty" if it contains only "dense" and/or
571 * "empty" columns/rows.
573 * stats [1]: (same as stats [0])
575 * stats [2]: number of garbage collections performed.
577 * stats [3]: status code. < 0 is an error code.
578 * > 1 is a warning or notice.
580 * 0 OK. Each column of the input matrix contained
581 * row indices in increasing order, with no
582 * duplicates.
584 * 1 OK, but columns of input matrix were jumbled
585 * (unsorted columns or duplicate entries). Symamd
586 * had to do some extra work to sort the matrix
587 * first and remove duplicate entries, but it
588 * still was able to return a valid permutation
589 * (return value of symamd was TRUE).
591 * stats [4]: highest numbered column that
592 * is unsorted or has duplicate
593 * entries.
594 * stats [5]: last seen duplicate or
595 * unsorted row index.
596 * stats [6]: number of duplicate or
597 * unsorted row indices.
599 * -1 A is a null pointer
601 * -2 p is a null pointer
603 * -3 (unused, see colamd.c)
605 * -4 n is negative
607 * stats [4]: n
609 * -5 number of nonzeros in matrix is negative
611 * stats [4]: # of nonzeros (p [n]).
613 * -6 p [0] is nonzero
615 * stats [4]: p [0]
617 * -7 (unused)
619 * -8 a column has a negative number of entries
621 * stats [4]: column with < 0 entries
622 * stats [5]: number of entries in col
624 * -9 a row index is out of bounds
626 * stats [4]: column with bad row index
627 * stats [5]: bad row index
628 * stats [6]: n_row, # of rows of matrx
630 * -10 out of memory (unable to allocate temporary
631 * workspace for M or count arrays using the
632 * "allocate" routine passed into symamd).
634 * Future versions may return more statistics in the stats array.
635 * @param allocate pointer to calloc
636 * @param release pointer to free
637 * @return TRUE if OK, FALSE otherwise
639 public static int symamd (int n, int[] A, int[] p, int[] perm,
640 double[] knobs, int[] stats)
642 /* === Local variables ============================================== */
644 int[] count ; /* length of each column of M, and col pointer*/
645 int[] mark ; /* mark array for finding duplicate entries */
646 int[] M ; /* row indices of matrix M */
647 int Mlen ; /* length of M */
648 int n_row ; /* number of rows in M */
649 int nnz ; /* number of entries in A */
650 int i ; /* row index of A */
651 int j ; /* column index of A */
652 int k ; /* row index of M */
653 int mnz ; /* number of nonzeros in M */
654 int pp ; /* index into a column of A */
655 int last_row ; /* last row seen in the current column */
656 int length ; /* number of nonzeros in a column */
658 double[] cknobs = new double[COLAMD_KNOBS] ; /* knobs for colamd */
659 double[] default_knobs = new double[COLAMD_KNOBS] ; /* default knobs for colamd */
661 if (!NDEBUG)
663 // colamd_get_debug ("symamd") ;
666 /* === Check the input arguments ==================================== */
668 if (stats == null)
670 DEBUG0 ("symamd: stats not present\n") ;
671 return (FALSE) ;
673 for (i = 0 ; i < COLAMD_STATS ; i++)
675 stats [i] = 0 ;
677 stats [COLAMD_STATUS] = COLAMD_OK ;
678 stats [COLAMD_INFO1] = -1 ;
679 stats [COLAMD_INFO2] = -1 ;
681 if (A == null)
683 stats [COLAMD_STATUS] = COLAMD_ERROR_A_not_present ;
684 DEBUG0 ("symamd: A not present\n") ;
685 return (FALSE) ;
688 if (p == null) /* p is not present */
690 stats [COLAMD_STATUS] = COLAMD_ERROR_p_not_present ;
691 DEBUG0 ("symamd: p not present\n") ;
692 return (FALSE) ;
695 if (n < 0) /* n must be >= 0 */
697 stats [COLAMD_STATUS] = COLAMD_ERROR_ncol_negative ;
698 stats [COLAMD_INFO1] = n ;
699 DEBUG0 ("symamd: n negative %d\n", n) ;
700 return (FALSE) ;
703 nnz = p [n] ;
704 if (nnz < 0) /* nnz must be >= 0 */
706 stats [COLAMD_STATUS] = COLAMD_ERROR_nnz_negative ;
707 stats [COLAMD_INFO1] = nnz ;
708 DEBUG0 ("symamd: number of entries negative %d\n", nnz) ;
709 return (FALSE) ;
712 if (p [0] != 0)
714 stats [COLAMD_STATUS] = COLAMD_ERROR_p0_nonzero ;
715 stats [COLAMD_INFO1] = p [0] ;
716 DEBUG0 ("symamd: p[0] not zero %d\n", p [0]) ;
717 return (FALSE) ;
720 /* === If no knobs, set default knobs =============================== */
722 if (knobs == null)
724 COLAMD_set_defaults (default_knobs) ;
725 knobs = default_knobs ;
728 /* === Allocate count and mark ====================================== */
732 count = new int [n+1] ;
733 } catch (OutOfMemoryError e) {
734 stats [COLAMD_STATUS] = COLAMD_ERROR_out_of_memory ;
735 DEBUG0 ("symamd: allocate count (size %d) failed\n", n+1) ;
736 return (FALSE) ;
741 mark = new int[n+1] ;
742 } catch (OutOfMemoryError e) {
743 stats [COLAMD_STATUS] = COLAMD_ERROR_out_of_memory ;
744 count = null ;
745 DEBUG0 ("symamd: allocate mark (size %d) failed\n", n+1) ;
746 return (FALSE) ;
749 /* === Compute column counts of M, check if A is valid ============== */
751 stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
753 for (i = 0 ; i < n ; i++)
755 mark [i] = -1 ;
758 for (j = 0 ; j < n ; j++)
760 last_row = -1 ;
762 length = p [j+1] - p [j] ;
763 if (length < 0)
765 /* column pointers must be non-decreasing */
766 stats [COLAMD_STATUS] = COLAMD_ERROR_col_length_negative ;
767 stats [COLAMD_INFO1] = j ;
768 stats [COLAMD_INFO2] = length ;
769 count = null ;
770 mark = null ;
771 DEBUG0 ("symamd: col %d negative length %d\n", j, length) ;
772 return (FALSE) ;
775 for (pp = p [j] ; pp < p [j+1] ; pp++)
777 i = A [pp] ;
778 if (i < 0 || i >= n)
780 /* row index i, in column j, is out of bounds */
781 stats [COLAMD_STATUS] = COLAMD_ERROR_row_index_out_of_bounds ;
782 stats [COLAMD_INFO1] = j ;
783 stats [COLAMD_INFO2] = i ;
784 stats [COLAMD_INFO3] = n ;
785 count = null ;
786 mark = null ;
787 DEBUG0 ("symamd: row %d col %d out of bounds\n", i, j) ;
788 return (FALSE) ;
791 if (i <= last_row || mark [i] == j)
793 /* row index is unsorted or repeated (or both), thus col */
794 /* is jumbled. This is a notice, not an error condition. */
795 stats [COLAMD_STATUS] = COLAMD_OK_BUT_JUMBLED ;
796 stats [COLAMD_INFO1] = j ;
797 stats [COLAMD_INFO2] = i ;
798 (stats [COLAMD_INFO3]) ++ ;
799 DEBUG1 ("symamd: row %d col %d unsorted/duplicate\n", i, j) ;
802 if (i > j && mark [i] != j)
804 /* row k of M will contain column indices i and j */
805 count [i]++ ;
806 count [j]++ ;
809 /* mark the row as having been seen in this column */
810 mark [i] = j ;
812 last_row = i ;
816 /* === Compute column pointers of M ================================= */
818 /* use output permutation, perm, for column pointers of M */
819 perm [0] = 0 ;
820 for (j = 1 ; j <= n ; j++)
822 perm [j] = perm [j-1] + count [j-1] ;
824 for (j = 0 ; j < n ; j++)
826 count [j] = perm [j] ;
829 /* === Construct M ================================================== */
831 mnz = perm [n] ;
832 n_row = mnz / 2 ;
833 Mlen = COLAMD_recommended (mnz, n_row, n) ;
836 M = new int [Mlen] ;
837 DEBUG0 ("symamd: M is %d-by-%d with %d entries, Mlen = %d\n",
838 n_row, n, mnz, Mlen) ;
839 } catch (OutOfMemoryError e) {
840 stats [COLAMD_STATUS] = COLAMD_ERROR_out_of_memory ;
841 count = null ;
842 mark = null;
843 DEBUG0 ("symamd: allocate M (size %g) failed\n", (double) Mlen) ;
844 return (FALSE) ;
847 k = 0 ;
849 if (stats [COLAMD_STATUS] == COLAMD_OK)
851 /* Matrix is OK */
852 for (j = 0 ; j < n ; j++)
854 ASSERT (p [j+1] - p [j] >= 0) ;
855 for (pp = p [j] ; pp < p [j+1] ; pp++)
857 i = A [pp] ;
858 ASSERT (i >= 0 && i < n) ;
859 if (i > j)
861 /* row k of M contains column indices i and j */
862 M [count [i]++] = k ;
863 M [count [j]++] = k ;
864 k++ ;
869 else
871 /* Matrix is jumbled. Do not add duplicates to M. Unsorted cols OK. */
872 DEBUG0 ("symamd: Duplicates in A.\n") ;
873 for (i = 0 ; i < n ; i++)
875 mark [i] = -1 ;
877 for (j = 0 ; j < n ; j++)
879 ASSERT (p [j+1] - p [j] >= 0) ;
880 for (pp = p [j] ; pp < p [j+1] ; pp++)
882 i = A [pp] ;
883 ASSERT (i >= 0 && i < n) ;
884 if (i > j && mark [i] != j)
886 /* row k of M contains column indices i and j */
887 M [count [i]++] = k ;
888 M [count [j]++] = k ;
889 k++ ;
890 mark [i] = j ;
896 /* count and mark no longer needed */
897 count = null ;
898 mark = null ;
899 ASSERT (k == n_row) ;
901 /* === Adjust the knobs for M ======================================= */
903 for (i = 0 ; i < COLAMD_KNOBS ; i++)
905 cknobs [i] = knobs [i] ;
908 /* there are no dense rows in M */
909 cknobs [COLAMD_DENSE_ROW] = -1 ;
910 cknobs [COLAMD_DENSE_COL] = knobs [COLAMD_DENSE_ROW] ;
912 /* === Order the columns of M ======================================= */
914 colamd (n_row, n, Mlen, M, perm, cknobs, stats) ;
916 /* Note that the output permutation is now in perm */
918 /* === get the statistics for symamd from colamd ==================== */
920 /* a dense column in colamd means a dense row and col in symamd */
921 stats [COLAMD_DENSE_ROW] = stats [COLAMD_DENSE_COL] ;
923 /* === Free M ======================================================= */
925 M = null ;
926 DEBUG0 ("symamd: done.\n") ;
927 return (TRUE) ;
930 /* ========================================================================== */
931 /* === colamd =============================================================== */
932 /* ========================================================================== */
935 The colamd routine computes a column ordering Q of a sparse matrix
936 A such that the LU factorization P(AQ) = LU remains sparse, where P is
937 selected via partial pivoting. The routine can also be viewed as
938 providing a permutation Q such that the Cholesky factorization
939 (AQ)'(AQ) = LL' remains sparse.
943 * Computes a column ordering (Q) of A such that P(AQ)=LU or
944 * (AQ)'AQ=LL' have less fill-in and require fewer floating point
945 * operations than factorizing the unpermuted matrix A or A'A,
946 * respectively.
948 * @param n_row number of rows in A. Restriction: n_row >= 0.
949 * Colamd returns FALSE if n_row is negative.
950 * @param n_col number of columns in A. Restriction: n_col >= 0.
951 * Colamd returns FALSE if n_col is negative.
952 * @param Alen length of A. Restriction (see note):
953 * Alen >= 2*nnz + 6*(n_col+1) + 4*(n_row+1) + n_col
954 * Colamd returns FALSE if these conditions are not met.
956 * Note: this restriction makes an modest assumption regarding
957 * the size of the two typedef's structures in colamd.h.
958 * We do, however, guarantee that
960 * Alen >= colamd_recommended (nnz, n_row, n_col)
962 * will be sufficient. Note: the macro version does not check
963 * for integer overflow, and thus is not recommended. Use
964 * the colamd_recommended routine instead.
965 * @param A row indices of A.
967 * A is an integer array of size Alen. Alen must be at least as
968 * large as the bare minimum value given above, but this is very
969 * low, and can result in excessive run time. For best
970 * performance, we recommend that Alen be greater than or equal to
971 * colamd_recommended (nnz, n_row, n_col), which adds
972 * nnz/5 to the bare minimum value given above.
974 * On input, the row indices of the entries in column c of the
975 * matrix are held in A [(p [c]) ... (p [c+1]-1)]. The row indices
976 * in a given column c need not be in ascending order, and
977 * duplicate row indices may be be present. However, colamd will
978 * work a little faster if both of these conditions are met
979 * (Colamd puts the matrix into this format, if it finds that the
980 * the conditions are not met).
982 * The matrix is 0-based. That is, rows are in the range 0 to
983 * n_row-1, and columns are in the range 0 to n_col-1. Colamd
984 * returns FALSE if any row index is out of range.
986 * The contents of A are modified during ordering, and are
987 * undefined on output.
988 * @param p pointers to columns in A.
990 * p is an integer array of size n_col+1. On input, it holds the
991 * "pointers" for the column form of the matrix A. Column c of
992 * the matrix A is held in A [(p [c]) ... (p [c+1]-1)]. The first
993 * entry, p [0], must be zero, and p [c] <= p [c+1] must hold
994 * for all c in the range 0 to n_col-1. The value p [n_col] is
995 * thus the total number of entries in the pattern of the matrix A.
996 * Colamd returns FALSE if these conditions are not met.
998 * On output, if colamd returns TRUE, the array p holds the column
999 * permutation (Q, for P(AQ)=LU or (AQ)'(AQ)=LL'), where p [0] is
1000 * the first column index in the new ordering, and p [n_col-1] is
1001 * the last. That is, p [k] = j means that column j of A is the
1002 * kth pivot column, in AQ, where k is in the range 0 to n_col-1
1003 * (p [0] = j means that column j of A is the first column in AQ).
1005 * If colamd returns FALSE, then no permutation is returned, and
1006 * p is undefined on output.
1007 * @param knobs parameters (uses defaults if NULL)
1008 * @param stats output statistics and error codes.
1010 * Statistics on the ordering, and error status.
1011 * Colamd returns FALSE if stats is not present.
1013 * stats [0]: number of dense or empty rows ignored.
1015 * stats [1]: number of dense or empty columns ignored (and
1016 * ordered last in the output permutation p)
1017 * Note that a row can become "empty" if it
1018 * contains only "dense" and/or "empty" columns,
1019 * and similarly a column can become "empty" if it
1020 * only contains "dense" and/or "empty" rows.
1022 * stats [2]: number of garbage collections performed.
1023 * This can be excessively high if Alen is close
1024 * to the minimum required value.
1026 * stats [3]: status code. < 0 is an error code.
1027 * > 1 is a warning or notice.
1029 * 0 OK. Each column of the input matrix contained
1030 * row indices in increasing order, with no
1031 * duplicates.
1033 * 1 OK, but columns of input matrix were jumbled
1034 * (unsorted columns or duplicate entries). Colamd
1035 * had to do some extra work to sort the matrix
1036 * first and remove duplicate entries, but it
1037 * still was able to return a valid permutation
1038 * (return value of colamd was TRUE).
1040 * stats [4]: highest numbered column that
1041 * is unsorted or has duplicate
1042 * entries.
1043 * stats [5]: last seen duplicate or
1044 * unsorted row index.
1045 * stats [6]: number of duplicate or
1046 * unsorted row indices.
1048 * -1 A is a null pointer
1050 * -2 p is a null pointer
1052 * -3 n_row is negative
1054 * stats [4]: n_row
1056 * -4 n_col is negative
1058 * stats [4]: n_col
1060 * -5 number of nonzeros in matrix is negative
1062 * stats [4]: number of nonzeros, p [n_col]
1064 * -6 p [0] is nonzero
1066 * stats [4]: p [0]
1068 * -7 A is too small
1070 * stats [4]: required size
1071 * stats [5]: actual size (Alen)
1073 * -8 a column has a negative number of entries
1075 * stats [4]: column with < 0 entries
1076 * stats [5]: number of entries in col
1078 * -9 a row index is out of bounds
1080 * stats [4]: column with bad row index
1081 * stats [5]: bad row index
1082 * stats [6]: n_row, # of rows of matrx
1084 * -10 (unused; see symamd.c)
1086 * -999 (unused; see symamd.c)
1088 * Future versions may return more statistics in the stats array.
1089 * @return TRUE if successful, FALSE otherwise
1091 public static int colamd (int n_row, int n_col, int Alen, int[] A,
1092 int[] p, double[] knobs, int[] stats)
1094 /* === Local variables ============================================== */
1096 int i ; /* loop index */
1097 int nnz ; /* nonzeros in A */
1098 int Row_size ; /* size of Row [], in integers */
1099 int Col_size ; /* size of Col [], in integers */
1100 int need ; /* minimum required length of A */
1101 Colamd_Row[] Row ; /* pointer into A of Row [0..n_row] array */
1102 Colamd_Col[] Col ; /* pointer into A of Col [0..n_col] array */
1103 int[] n_col2 = new int [1] ; /* number of non-dense, non-empty columns */
1104 int[] n_row2 = new int [1] ; /* number of non-dense, non-empty rows */
1105 int ngarbage ; /* number of garbage collections performed */
1106 int[] max_deg = new int [1] ; /* maximum row degree */
1107 double[] default_knobs = new double[COLAMD_KNOBS] ; /* default knobs array */
1108 int aggressive ; /* do aggressive absorption */
1109 int[] ok ;
1111 if (!NDEBUG)
1113 // colamd_get_debug ("colamd") ;
1116 /* === Check the input arguments ==================================== */
1118 if (stats == null)
1120 DEBUG0 ("colamd: stats not present\n") ;
1121 return (FALSE) ;
1123 for (i = 0 ; i < COLAMD_STATS ; i++)
1125 stats [i] = 0 ;
1127 stats [COLAMD_STATUS] = COLAMD_OK ;
1128 stats [COLAMD_INFO1] = -1 ;
1129 stats [COLAMD_INFO2] = -1 ;
1131 if (A == null) /* A is not present */
1133 stats [COLAMD_STATUS] = COLAMD_ERROR_A_not_present ;
1134 DEBUG0 ("colamd: A not present\n") ;
1135 return (FALSE) ;
1138 if (p == null) /* p is not present */
1140 stats [COLAMD_STATUS] = COLAMD_ERROR_p_not_present ;
1141 DEBUG0 ("colamd: p not present\n") ;
1142 return (FALSE) ;
1145 if (n_row < 0) /* n_row must be >= 0 */
1147 stats [COLAMD_STATUS] = COLAMD_ERROR_nrow_negative ;
1148 stats [COLAMD_INFO1] = n_row ;
1149 DEBUG0 ("colamd: nrow negative %d\n", n_row) ;
1150 return (FALSE) ;
1153 if (n_col < 0) /* n_col must be >= 0 */
1155 stats [COLAMD_STATUS] = COLAMD_ERROR_ncol_negative ;
1156 stats [COLAMD_INFO1] = n_col ;
1157 DEBUG0 ("colamd: ncol negative %d\n", n_col) ;
1158 return (FALSE) ;
1161 nnz = p [n_col] ;
1162 if (nnz < 0) /* nnz must be >= 0 */
1164 stats [COLAMD_STATUS] = COLAMD_ERROR_nnz_negative ;
1165 stats [COLAMD_INFO1] = nnz ;
1166 DEBUG0 ("colamd: number of entries negative %d\n", nnz) ;
1167 return (FALSE) ;
1170 if (p [0] != 0)
1172 stats [COLAMD_STATUS] = COLAMD_ERROR_p0_nonzero ;
1173 stats [COLAMD_INFO1] = p [0] ;
1174 DEBUG0 ("colamd: p[0] not zero %d\n", p [0]) ;
1175 return (FALSE) ;
1178 /* === If no knobs, set default knobs =============================== */
1180 if (knobs == null)
1182 COLAMD_set_defaults (default_knobs) ;
1183 knobs = default_knobs ;
1186 aggressive = (knobs [COLAMD_AGGRESSIVE] != FALSE) ? 1 : 0;
1188 /* === Allocate the Row and Col arrays from array A ================= */
1190 ok = new int [] { TRUE } ;
1191 Col_size = COLAMD_C (n_col, ok) ; /* size of Col array of structs */
1192 Row_size = COLAMD_R (n_row, ok) ; /* size of Row array of structs */
1194 /* need = 2*nnz + n_col + Col_size + Row_size ; */
1195 need = t_mult (nnz, 2, ok) ;
1196 need = t_add (need, n_col, ok) ;
1197 // need = t_add (need, Col_size, ok) ;
1198 // need = t_add (need, Row_size, ok) ;
1200 if ((ok[0] == 0) || need > (int) Alen || need > Int_MAX)
1202 /* not enough space in array A to perform the ordering */
1203 stats [COLAMD_STATUS] = COLAMD_ERROR_A_too_small ;
1204 stats [COLAMD_INFO1] = need ;
1205 stats [COLAMD_INFO2] = Alen ;
1206 DEBUG0 ("colamd: Need Alen >= %d, given only Alen = %d\n", need,Alen);
1207 return (FALSE) ;
1210 // Alen -= Col_size + Row_size ;
1211 Col = new Colamd_Col [Col_size] ; //A [Alen] ;
1212 Row = new Colamd_Row [Row_size] ; //A [Alen + Col_size] ;
1214 for (i = 0; i < Col_size; i++)
1215 Col [i] = new Colamd_Col() ;
1216 for (i = 0; i < Row_size; i++)
1217 Row [i] = new Colamd_Row() ;
1219 /* === Construct the row and column data structures ================= */
1221 if (init_rows_cols (n_row, n_col, Row, Col, A, p, stats) == 0)
1223 /* input matrix is invalid */
1224 DEBUG0 ("colamd: Matrix invalid\n") ;
1225 return (FALSE) ;
1228 /* === Initialize scores, kill dense rows/columns =================== */
1230 init_scoring (n_row, n_col, Row, Col, A, p, knobs,
1231 n_row2, n_col2, max_deg) ;
1233 /* === Order the supercolumns ======================================= */
1235 ngarbage = find_ordering (n_row, n_col, Alen, Row, Col, A, p,
1236 n_col2[0], max_deg[0], 2*nnz, aggressive) ;
1238 /* === Order the non-principal columns ============================== */
1240 order_children (n_col, Col, p) ;
1242 /* === Return statistics in stats =================================== */
1244 stats [COLAMD_DENSE_ROW] = n_row - n_row2[0] ;
1245 stats [COLAMD_DENSE_COL] = n_col - n_col2[0] ;
1246 stats [COLAMD_DEFRAG_COUNT] = ngarbage ;
1247 DEBUG0 ("colamd: done.\n") ;
1248 return (TRUE) ;
1252 /* ========================================================================== */
1253 /* === colamd_report ======================================================== */
1254 /* ========================================================================== */
1256 public static void COLAMD_report (int[] stats)
1258 print_report ("colamd", stats) ;
1262 /* ========================================================================== */
1263 /* === symamd_report ======================================================== */
1264 /* ========================================================================== */
1266 public static void SYMAMD_report (int[] stats)
1268 print_report ("symamd", stats) ;
1273 /* ========================================================================== */
1274 /* === NON-USER-CALLABLE ROUTINES: ========================================== */
1275 /* ========================================================================== */
1277 /* There are no user-callable routines beyond this point in the file */
1280 /* ========================================================================== */
1281 /* === init_rows_cols ======================================================= */
1282 /* ========================================================================== */
1285 * Takes the column form of the matrix in A and creates the row form of the
1286 * matrix. Also, row and column attributes are stored in the Col and Row
1287 * structs. If the columns are un-sorted or contain duplicate row indices,
1288 * this routine will also sort and remove duplicate row indices from the
1289 * column form of the matrix. Returns FALSE if the matrix is invalid,
1290 * TRUE otherwise. Not user-callable.
1292 * @param n_row number of rows of A
1293 * @param n_col number of columns of A
1294 * @param Row of size n_row+1
1295 * @param Col of size n_col+1
1296 * @param A row indices of A, of size Alen
1297 * @param p pointers to columns in A, of size n_col+1
1298 * @param stats colamd statistics
1299 * @return TRUE if OK, or FALSE otherwise
1301 private static int init_rows_cols (int n_row, int n_col, Colamd_Row[] Row,
1302 Colamd_Col[] Col, int[] A, int[] p, int[] stats)
1304 /* === Local variables ============================================== */
1306 int i ;
1307 int col ; /* a column index */
1308 int row ; /* a row index */
1309 int cp ; /* a column pointer */
1310 int cp_end ; /* a pointer to the end of a column */
1311 int rp ; /* a row pointer */
1312 int rp_end ; /* a pointer to the end of a row */
1313 int last_row ; /* previous row */
1315 /* === Initialize columns, and check column pointers ================ */
1317 for (col = 0 ; col < n_col ; col++)
1319 Col [col].start = p [col] ;
1320 Col [col].length = p [col+1] - p [col] ;
1322 if (Col [col].length < 0)
1324 /* column pointers must be non-decreasing */
1325 stats [COLAMD_STATUS] = COLAMD_ERROR_col_length_negative ;
1326 stats [COLAMD_INFO1] = col ;
1327 stats [COLAMD_INFO2] = Col [col].length ;
1328 DEBUG0 ("colamd: col %d length %d < 0\n", col, Col [col].length) ;
1329 return (FALSE) ;
1332 Col [col].thickness = 1 ;
1333 Col [col].score = 0 ;
1334 Col [col].prev = EMPTY ;
1335 Col [col].degree_next = EMPTY ;
1338 /* p [0..n_col] no longer needed, used as "head" in subsequent routines */
1340 /* === Scan columns, compute row degrees, and check row indices ===== */
1342 stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
1344 for (row = 0 ; row < n_row ; row++)
1346 Row [row].length = 0 ;
1347 Row [row].mark = -1 ;
1350 for (col = 0 ; col < n_col ; col++)
1352 last_row = -1 ;
1354 cp = p [col] ;
1355 cp_end = p [col+1] ;
1357 while (cp < cp_end)
1359 row = A [cp++] ;
1361 /* make sure row indices within range */
1362 if (row < 0 || row >= n_row)
1364 stats [COLAMD_STATUS] = COLAMD_ERROR_row_index_out_of_bounds ;
1365 stats [COLAMD_INFO1] = col ;
1366 stats [COLAMD_INFO2] = row ;
1367 stats [COLAMD_INFO3] = n_row ;
1368 DEBUG0 ("colamd: row %d col %d out of bounds\n", row, col) ;
1369 return (FALSE) ;
1372 if (row <= last_row || Row [row].mark == col)
1374 /* row index are unsorted or repeated (or both), thus col */
1375 /* is jumbled. This is a notice, not an error condition. */
1376 stats [COLAMD_STATUS] = COLAMD_OK_BUT_JUMBLED ;
1377 stats [COLAMD_INFO1] = col ;
1378 stats [COLAMD_INFO2] = row ;
1379 (stats [COLAMD_INFO3]) ++ ;
1380 DEBUG1 ("colamd: row %d col %d unsorted/duplicate\n",row,col);
1383 if (Row [row].mark != col)
1385 Row [row].length++ ;
1387 else
1389 /* this is a repeated entry in the column, */
1390 /* it will be removed */
1391 Col [col].length-- ;
1394 /* mark the row as having been seen in this column */
1395 Row [row].mark = col ;
1397 last_row = row ;
1401 /* === Compute row pointers ========================================= */
1403 /* row form of the matrix starts directly after the column */
1404 /* form of matrix in A */
1405 Row [0].start = p [n_col] ;
1406 Row [0].p = Row [0].start ;
1407 Row [0].mark = -1 ;
1408 for (row = 1 ; row < n_row ; row++)
1410 Row [row].start = Row [row-1].start + Row [row-1].length ;
1411 Row [row].p = Row [row].start ;
1412 Row [row].mark = -1 ;
1415 /* === Create row form ============================================== */
1417 if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1419 /* if cols jumbled, watch for repeated row indices */
1420 for (col = 0 ; col < n_col ; col++)
1422 cp = p [col] ;
1423 cp_end = p [col+1] ;
1425 while (cp < cp_end)
1427 row = A [cp++] ;
1429 if (Row [row].mark != col)
1431 A [(Row [row].p)++] = col ;
1432 Row [row].mark = col ;
1437 else
1439 /* if cols not jumbled, we don't need the mark (this is faster) */
1440 for (col = 0 ; col < n_col ; col++)
1442 cp = p [col] ;
1443 cp_end = p [col+1] ;
1444 while (cp < cp_end)
1446 A [(Row [A [cp++]].p)++] = col ;
1451 /* === Clear the row marks and set row degrees ====================== */
1453 for (row = 0 ; row < n_row ; row++)
1455 Row [row].mark = 0 ;
1456 Row [row].degree = Row [row].length ;
1459 /* === See if we need to re-create columns ========================== */
1461 if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1463 DEBUG0 ("colamd: reconstructing column form, matrix jumbled\n") ;
1465 if (!NDEBUG)
1467 /* make sure column lengths are correct */
1468 for (col = 0 ; col < n_col ; col++)
1470 p [col] = Col [col].length ;
1472 for (row = 0 ; row < n_row ; row++)
1474 rp = Row [row].start ;
1475 rp_end = rp + Row [row].length ;
1476 while (rp < rp_end)
1478 p [A [rp++]]-- ;
1481 for (col = 0 ; col < n_col ; col++)
1483 ASSERT (p [col] == 0) ;
1485 /* now p is all zero (different than when debugging is turned off) */
1486 } /* NDEBUG */
1488 /* === Compute col pointers ========================================= */
1490 /* col form of the matrix starts at A [0]. */
1491 /* Note, we may have a gap between the col form and the row */
1492 /* form if there were duplicate entries, if so, it will be */
1493 /* removed upon the first garbage collection */
1494 Col [0].start = 0 ;
1495 p [0] = Col [0].start ;
1496 for (col = 1 ; col < n_col ; col++)
1498 /* note that the lengths here are for pruned columns, i.e. */
1499 /* no duplicate row indices will exist for these columns */
1500 Col [col].start = Col [col-1].start + Col [col-1].length ;
1501 p [col] = Col [col].start ;
1504 /* === Re-create col form =========================================== */
1506 for (row = 0 ; row < n_row ; row++)
1508 rp = Row [row].start ;
1509 rp_end = rp + Row [row].length ;
1510 while (rp < rp_end)
1512 A [(p [A [rp++]])++] = row ;
1517 /* === Done. Matrix is not (or no longer) jumbled ================== */
1519 return (TRUE) ;
1523 /* ========================================================================== */
1524 /* === init_scoring ========================================================= */
1525 /* ========================================================================== */
1528 * Kills dense or empty columns and rows, calculates an initial score for
1529 * each column, and places all columns in the degree lists.init_rows_cols
1531 * @param n_row number of rows of A
1532 * @param n_col number of columns of A
1533 * @param Row of size n_row+1
1534 * @param Col of size n_col+1
1535 * @param A column form and row form of A
1536 * @param head of size n_col+1
1537 * @param knobs parameters
1538 * @param p_n_row2 size 1, number of non-dense, non-empty rows
1539 * @param p_n_col2 size 1, number of non-dense, non-empty columns
1540 * @param p_max_deg size 1, maximum row degree
1542 private static void init_scoring (int n_row, int n_col, Colamd_Row[] Row,
1543 Colamd_Col[] Col, int[] A, int[] head, double[] knobs,
1544 int[] p_n_row2, int[] p_n_col2, int[] p_max_deg)
1546 /* === Local variables ============================================== */
1548 int c ; /* a column index */
1549 int r, row ; /* a row index */
1550 int cp ; /* a column pointer */
1551 int deg ; /* degree of a row or column */
1552 int cp_end ; /* a pointer to the end of a column */
1553 int new_cp ; /* new column pointer */
1554 int col_length ; /* length of pruned column */
1555 int score ; /* current column score */
1556 int n_col2 ; /* number of non-dense, non-empty columns */
1557 int n_row2 ; /* number of non-dense, non-empty rows */
1558 int dense_row_count ; /* remove rows with more entries than this */
1559 int dense_col_count ; /* remove cols with more entries than this */
1560 int min_score ; /* smallest column score */
1561 int max_deg ; /* maximum row degree */
1562 int next_col ; /* Used to add to degree list.*/
1564 int debug_count = 0 ; /* debug only. */
1566 /* === Extract knobs ================================================ */
1568 /* Note: if knobs contains a NaN, this is undefined: */
1569 if (knobs [COLAMD_DENSE_ROW] < 0)
1571 /* only remove completely dense rows */
1572 dense_row_count = n_col-1 ;
1574 else
1576 dense_row_count = DENSE_DEGREE (knobs [COLAMD_DENSE_ROW], n_col) ;
1578 if (knobs [COLAMD_DENSE_COL] < 0)
1580 /* only remove completely dense columns */
1581 dense_col_count = n_row-1 ;
1583 else
1585 dense_col_count = DENSE_DEGREE (knobs [COLAMD_DENSE_COL],
1586 MIN (n_row, n_col)) ;
1589 DEBUG1 ("colamd: densecount: %d %d\n", dense_row_count, dense_col_count) ;
1590 max_deg = 0 ;
1591 n_col2 = n_col ;
1592 n_row2 = n_row ;
1594 /* === Kill empty columns =========================================== */
1596 /* Put the empty columns at the end in their natural order, so that LU */
1597 /* factorization can proceed as far as possible. */
1598 for (c = n_col-1 ; c >= 0 ; c--)
1600 deg = Col [c].length ;
1601 if (deg == 0)
1603 /* this is a empty column, kill and order it last */
1604 Col [c].order = --n_col2 ;
1605 KILL_PRINCIPAL_COL (Col, c) ;
1608 DEBUG1 ("colamd: null columns killed: %d\n", n_col - n_col2) ;
1610 /* === Kill dense columns =========================================== */
1612 /* Put the dense columns at the end, in their natural order */
1613 for (c = n_col-1 ; c >= 0 ; c--)
1615 /* skip any dead columns */
1616 if (COL_IS_DEAD (Col, c))
1618 continue ;
1620 deg = Col [c].length ;
1621 if (deg > dense_col_count)
1623 /* this is a dense column, kill and order it last */
1624 Col [c].order = --n_col2 ;
1625 /* decrement the row degrees */
1626 cp = Col [c].start ;
1627 cp_end = cp + Col [c].length ;
1628 while (cp < cp_end)
1630 Row [A [cp++]].degree-- ;
1632 KILL_PRINCIPAL_COL (Col, c) ;
1635 DEBUG1 ("colamd: Dense and null columns killed: %d\n", n_col - n_col2) ;
1637 /* === Kill dense and empty rows ==================================== */
1639 for (r = 0 ; r < n_row ; r++)
1641 deg = Row [r].degree ;
1642 ASSERT (deg >= 0 && deg <= n_col) ;
1643 if (deg > dense_row_count || deg == 0)
1645 /* kill a dense or empty row */
1646 KILL_ROW (Row, r) ;
1647 --n_row2 ;
1649 else
1651 /* keep track of max degree of remaining rows */
1652 max_deg = MAX (max_deg, deg) ;
1655 DEBUG1 ("colamd: Dense and null rows killed: %d\n", n_row - n_row2) ;
1657 /* === Compute initial column scores ================================ */
1659 /* At this point the row degrees are accurate. They reflect the number */
1660 /* of "live" (non-dense) columns in each row. No empty rows exist. */
1661 /* Some "live" columns may contain only dead rows, however. These are */
1662 /* pruned in the code below. */
1664 /* now find the initial matlab score for each column */
1665 for (c = n_col-1 ; c >= 0 ; c--)
1667 /* skip dead column */
1668 if (COL_IS_DEAD (Col, c))
1670 continue ;
1672 score = 0 ;
1673 cp = Col [c].start ;
1674 new_cp = cp ;
1675 cp_end = cp + Col [c].length ;
1676 while (cp < cp_end)
1678 /* get a row */
1679 row = A [cp++] ;
1680 /* skip if dead */
1681 if (ROW_IS_DEAD (Row, row))
1683 continue ;
1685 /* compact the column */
1686 A [new_cp++] = row ;
1687 /* add row's external degree */
1688 score += Row [row].degree - 1 ;
1689 /* guard against integer overflow */
1690 score = MIN (score, n_col) ;
1692 /* determine pruned column length */
1693 col_length = (new_cp - Col [c].start) ;
1694 if (col_length == 0)
1696 /* a newly-made null column (all rows in this col are "dense" */
1697 /* and have already been killed) */
1698 DEBUG2 ("Newly null killed: %d\n", c) ;
1699 Col [c].order = --n_col2 ;
1700 KILL_PRINCIPAL_COL (Col, c) ;
1702 else
1704 /* set column length and set score */
1705 ASSERT (score >= 0) ;
1706 ASSERT (score <= n_col) ;
1707 Col [c].length = col_length ;
1708 Col [c].score = score ;
1711 DEBUG1 ("colamd: Dense, null, and newly-null columns killed: %d\n",
1712 n_col-n_col2) ;
1714 /* At this point, all empty rows and columns are dead. All live columns */
1715 /* are "clean" (containing no dead rows) and simplicial (no supercolumns */
1716 /* yet). Rows may contain dead columns, but all live rows contain at */
1717 /* least one live column. */
1719 if (!NDEBUG)
1721 debug_structures (n_row, n_col, Row, Col, A, n_col2) ;
1724 /* === Initialize degree lists ========================================== */
1726 if (!NDEBUG)
1728 debug_count = 0 ;
1731 /* clear the hash buckets */
1732 for (c = 0 ; c <= n_col ; c++)
1734 head [c] = EMPTY ;
1736 min_score = n_col ;
1737 /* place in reverse order, so low column indices are at the front */
1738 /* of the lists. This is to encourage natural tie-breaking */
1739 for (c = n_col-1 ; c >= 0 ; c--)
1741 /* only add principal columns to degree lists */
1742 if (COL_IS_ALIVE (Col, c))
1744 DEBUG4 ("place %d score %d minscore %d ncol %d\n",
1745 c, Col [c].score, min_score, n_col) ;
1747 /* === Add columns score to DList =============================== */
1749 score = Col [c].score ;
1751 ASSERT (min_score >= 0) ;
1752 ASSERT (min_score <= n_col) ;
1753 ASSERT (score >= 0) ;
1754 ASSERT (score <= n_col) ;
1755 ASSERT (head [score] >= EMPTY) ;
1757 /* now add this column to dList at proper score location */
1758 next_col = head [score] ;
1759 Col [c].prev = EMPTY ;
1760 Col [c].degree_next = next_col ;
1762 /* if there already was a column with the same score, set its */
1763 /* previous pointer to this new column */
1764 if (next_col != EMPTY)
1766 Col [next_col].prev = c ;
1768 head [score] = c ;
1770 /* see if this score is less than current min */
1771 min_score = MIN (min_score, score) ;
1773 if (!NDEBUG)
1775 debug_count++ ;
1781 if (!NDEBUG)
1783 DEBUG1 ("colamd: Live cols %d out of %d, non-princ: %d\n",
1784 debug_count, n_col, n_col-debug_count) ;
1785 ASSERT (debug_count == n_col2) ;
1786 debug_deg_lists (n_row, n_col, Row, Col, head, min_score, n_col2, max_deg) ;
1787 } /* NDEBUG */
1789 /* === Return number of remaining columns, and max row degree ======= */
1791 p_n_col2[0] = n_col2 ;
1792 p_n_row2[0] = n_row2 ;
1793 p_max_deg[0] = max_deg ;
1797 /* ========================================================================== */
1798 /* === find_ordering ======================================================== */
1799 /* ========================================================================== */
1802 * Order the principal columns of the supercolumn form of the matrix
1803 * (no supercolumns on input). Uses a minimum approximate column minimum
1804 * degree ordering method. Not user-callable.
1806 * @param n_row number of rows of A
1807 * @param n_col number of columns of A
1808 * @param Alen size of A, 2*nnz + n_col or larger
1809 * @param Row of size n_row+1
1810 * @param Col of size n_col+1
1811 * @param A column form and row form of A
1812 * @param head of size n_col+1
1813 * @param n_col2 Remaining columns to order
1814 * @param max_deg Maximum row degree
1815 * @param pfree index of first free slot (2*nnz on entry)
1816 * @param aggressive
1817 * @return the number of garbage collections
1819 private static int find_ordering (int n_row, int n_col, int Alen,
1820 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int[] head, int n_col2,
1821 int max_deg, int pfree, int aggressive)
1823 /* === Local variables ============================================== */
1825 int k ; /* current pivot ordering step */
1826 int pivot_col ; /* current pivot column */
1827 int cp ; /* a column pointer */
1828 int rp ; /* a row pointer */
1829 int pivot_row ; /* current pivot row */
1830 int new_cp ; /* modified column pointer */
1831 int new_rp ; /* modified row pointer */
1832 int pivot_row_start ; /* pointer to start of pivot row */
1833 int pivot_row_degree ; /* number of columns in pivot row */
1834 int pivot_row_length ; /* number of supercolumns in pivot row */
1835 int pivot_col_score ; /* score of pivot column */
1836 int needed_memory ; /* free space needed for pivot row */
1837 int cp_end ; /* pointer to the end of a column */
1838 int rp_end ; /* pointer to the end of a row */
1839 int row ; /* a row index */
1840 int col ; /* a column index */
1841 int max_score ; /* maximum possible score */
1842 int cur_score ; /* score of current column */
1843 /* FIXME unsigned */ int hash ; /* hash value for supernode detection */
1844 int head_column ; /* head of hash bucket */
1845 int first_col ; /* first column in hash bucket */
1846 int tag_mark ; /* marker value for mark array */
1847 int row_mark ; /* Row [row].shared2.mark */
1848 int set_difference ; /* set difference size of row with pivot row */
1849 int min_score ; /* smallest column score */
1850 int col_thickness ; /* "thickness" (no. of columns in a supercol) */
1851 int max_mark ; /* maximum value of tag_mark */
1852 int pivot_col_thickness ; /* number of columns represented by pivot col */
1853 int prev_col ; /* Used by Dlist operations. */
1854 int next_col ; /* Used by Dlist operations. */
1855 int ngarbage ; /* number of garbage collections performed */
1857 int debug_d ; /* debug loop counter */
1858 int debug_step = 0 ; /* debug loop counter */
1860 /* === Initialization and clear mark ================================ */
1862 max_mark = Int_MAX - n_col ;
1863 tag_mark = clear_mark (0, max_mark, n_row, Row) ;
1864 min_score = 0 ;
1865 ngarbage = 0 ;
1866 DEBUG1 ("colamd: Ordering, n_col2=%d\n", n_col2) ;
1868 /* === Order the columns ============================================ */
1870 for (k = 0 ; k < n_col2 ; /* 'k' is incremented below */)
1873 if (!NDEBUG)
1875 if (debug_step % 100 == 0)
1877 DEBUG2 ("\n... Step k: %d out of n_col2: %d\n", k, n_col2) ;
1879 else
1881 DEBUG3 ("\n----------Step k: %d out of n_col2: %d\n", k, n_col2) ;
1883 debug_step++ ;
1884 debug_deg_lists (n_row, n_col, Row, Col, head,
1885 min_score, n_col2-k, max_deg) ;
1886 debug_matrix (n_row, n_col, Row, Col, A) ;
1887 } /* NDEBUG */
1889 /* === Select pivot column, and order it ============================ */
1891 /* make sure degree list isn't empty */
1892 ASSERT (min_score >= 0) ;
1893 ASSERT (min_score <= n_col) ;
1894 ASSERT (head [min_score] >= EMPTY) ;
1896 if (!NDEBUG)
1898 for (debug_d = 0 ; debug_d < min_score ; debug_d++)
1900 ASSERT (head [debug_d] == EMPTY) ;
1902 } /* NDEBUG */
1904 /* get pivot column from head of minimum degree list */
1905 while (head [min_score] == EMPTY && min_score < n_col)
1907 min_score++ ;
1909 pivot_col = head [min_score] ;
1910 ASSERT (pivot_col >= 0 && pivot_col <= n_col) ;
1911 next_col = Col [pivot_col].degree_next ;
1912 head [min_score] = next_col ;
1913 if (next_col != EMPTY)
1915 Col [next_col].prev = EMPTY ;
1918 ASSERT (COL_IS_ALIVE (Col, pivot_col)) ;
1920 /* remember score for defrag check */
1921 pivot_col_score = Col [pivot_col].score ;
1923 /* the pivot column is the kth column in the pivot order */
1924 Col [pivot_col].order = k ;
1926 /* increment order count by column thickness */
1927 pivot_col_thickness = Col [pivot_col].thickness ;
1928 k += pivot_col_thickness ;
1929 ASSERT (pivot_col_thickness > 0) ;
1930 DEBUG3 ("Pivot col: %d thick %d\n", pivot_col, pivot_col_thickness) ;
1932 /* === Garbage_collection, if necessary ============================= */
1934 needed_memory = MIN (pivot_col_score, n_col - k) ;
1935 if (pfree + needed_memory >= Alen)
1937 pfree = garbage_collection (n_row, n_col, Row, Col, A, pfree) ;
1938 ngarbage++ ;
1939 /* after garbage collection we will have enough */
1940 ASSERT (pfree + needed_memory < Alen) ;
1941 /* garbage collection has wiped out the Row[].shared2.mark array */
1942 tag_mark = clear_mark (0, max_mark, n_row, Row) ;
1944 if (!NDEBUG)
1946 debug_matrix (n_row, n_col, Row, Col, A) ;
1947 } /* NDEBUG */
1950 /* === Compute pivot row pattern ==================================== */
1952 /* get starting location for this new merged row */
1953 pivot_row_start = pfree ;
1955 /* initialize new row counts to zero */
1956 pivot_row_degree = 0 ;
1958 /* tag pivot column as having been visited so it isn't included */
1959 /* in merged pivot row */
1960 Col [pivot_col].thickness = -pivot_col_thickness ;
1962 /* pivot row is the union of all rows in the pivot column pattern */
1963 cp = Col [pivot_col].start ;
1964 cp_end = cp + Col [pivot_col].length ;
1965 while (cp < cp_end)
1967 /* get a row */
1968 row = A [cp++] ;
1969 DEBUG4 ("Pivot col pattern %d %d\n", ROW_IS_ALIVE (Row, row) ? 1 : 0, row) ;
1970 /* skip if row is dead */
1971 if (ROW_IS_ALIVE (Row, row))
1973 rp = Row [row].start ;
1974 rp_end = rp + Row [row].length ;
1975 while (rp < rp_end)
1977 /* get a column */
1978 col = A [rp++] ;
1979 /* add the column, if alive and untagged */
1980 col_thickness = Col [col].thickness ;
1981 if (col_thickness > 0 && COL_IS_ALIVE (Col, col))
1983 /* tag column in pivot row */
1984 Col [col].thickness = -col_thickness ;
1985 ASSERT (pfree < Alen) ;
1986 /* place column in pivot row */
1987 A [pfree++] = col ;
1988 pivot_row_degree += col_thickness ;
1994 /* clear tag on pivot column */
1995 Col [pivot_col].thickness = pivot_col_thickness ;
1996 max_deg = MAX (max_deg, pivot_row_degree) ;
1998 if (!NDEBUG)
2000 DEBUG3 ("check2\n") ;
2001 debug_mark (n_row, Row, tag_mark, max_mark) ;
2002 } /* NDEBUG */
2004 /* === Kill all rows used to construct pivot row ==================== */
2006 /* also kill pivot row, temporarily */
2007 cp = Col [pivot_col].start ;
2008 cp_end = cp + Col [pivot_col].length ;
2009 while (cp < cp_end)
2011 /* may be killing an already dead row */
2012 row = A [cp++] ;
2013 DEBUG3 ("Kill row in pivot col: %d\n", row) ;
2014 KILL_ROW (Row, row) ;
2017 /* === Select a row index to use as the new pivot row =============== */
2019 pivot_row_length = pfree - pivot_row_start ;
2020 if (pivot_row_length > 0)
2022 /* pick the "pivot" row arbitrarily (first row in col) */
2023 pivot_row = A [Col [pivot_col].start] ;
2024 DEBUG3 ("Pivotal row is %d\n", pivot_row) ;
2026 else
2028 /* there is no pivot row, since it is of zero length */
2029 pivot_row = EMPTY ;
2030 ASSERT (pivot_row_length == 0) ;
2032 ASSERT (Col [pivot_col].length > 0 || pivot_row_length == 0) ;
2034 /* === Approximate degree computation =============================== */
2036 /* Here begins the computation of the approximate degree. The column */
2037 /* score is the sum of the pivot row "length", plus the size of the */
2038 /* set differences of each row in the column minus the pattern of the */
2039 /* pivot row itself. The column ("thickness") itself is also */
2040 /* excluded from the column score (we thus use an approximate */
2041 /* external degree). */
2043 /* The time taken by the following code (compute set differences, and */
2044 /* add them up) is proportional to the size of the data structure */
2045 /* being scanned - that is, the sum of the sizes of each column in */
2046 /* the pivot row. Thus, the amortized time to compute a column score */
2047 /* is proportional to the size of that column (where size, in this */
2048 /* context, is the column "length", or the number of row indices */
2049 /* in that column). The number of row indices in a column is */
2050 /* monotonically non-decreasing, from the length of the original */
2051 /* column on input to colamd. */
2053 /* === Compute set differences ====================================== */
2055 DEBUG3 ("** Computing set differences phase. **\n") ;
2057 /* pivot row is currently dead - it will be revived later. */
2059 DEBUG3 ("Pivot row: ") ;
2060 /* for each column in pivot row */
2061 rp = pivot_row_start ;
2062 rp_end = rp + pivot_row_length ;
2063 while (rp < rp_end)
2065 col = A [rp++] ;
2066 ASSERT (COL_IS_ALIVE (Col, col) && col != pivot_col) ;
2067 DEBUG3 ("Col: %d\n", col) ;
2069 /* clear tags used to construct pivot row pattern */
2070 col_thickness = -Col [col].thickness ;
2071 ASSERT (col_thickness > 0) ;
2072 Col [col].thickness = col_thickness ;
2074 /* === Remove column from degree list =========================== */
2076 cur_score = Col [col].score ;
2077 prev_col = Col [col].prev ;
2078 next_col = Col [col].degree_next ;
2079 ASSERT (cur_score >= 0) ;
2080 ASSERT (cur_score <= n_col) ;
2081 ASSERT (cur_score >= EMPTY) ;
2082 if (prev_col == EMPTY)
2084 head [cur_score] = next_col ;
2086 else
2088 Col [prev_col].degree_next = next_col ;
2090 if (next_col != EMPTY)
2092 Col [next_col].prev = prev_col ;
2095 /* === Scan the column ========================================== */
2097 cp = Col [col].start ;
2098 cp_end = cp + Col [col].length ;
2099 while (cp < cp_end)
2101 /* get a row */
2102 row = A [cp++] ;
2103 row_mark = Row [row].mark ;
2104 /* skip if dead */
2105 if (ROW_IS_MARKED_DEAD (row_mark))
2107 continue ;
2109 ASSERT (row != pivot_row) ;
2110 set_difference = row_mark - tag_mark ;
2111 /* check if the row has been seen yet */
2112 if (set_difference < 0)
2114 ASSERT (Row [row].degree <= max_deg) ;
2115 set_difference = Row [row].degree ;
2117 /* subtract column thickness from this row's set difference */
2118 set_difference -= col_thickness ;
2119 ASSERT (set_difference >= 0) ;
2120 /* absorb this row if the set difference becomes zero */
2121 if (set_difference == 0 && aggressive != 0)
2123 DEBUG3 ("aggressive absorption. Row: %d\n", row) ;
2124 KILL_ROW (Row, row) ;
2126 else
2128 /* save the new mark */
2129 Row [row].mark = set_difference + tag_mark ;
2134 if (!NDEBUG)
2136 debug_deg_lists (n_row, n_col, Row, Col, head,
2137 min_score, n_col2-k-pivot_row_degree, max_deg) ;
2138 } /* NDEBUG */
2140 /* === Add up set differences for each column ======================= */
2142 DEBUG3 ("** Adding set differences phase. **\n") ;
2144 /* for each column in pivot row */
2145 rp = pivot_row_start ;
2146 rp_end = rp + pivot_row_length ;
2147 while (rp < rp_end)
2149 /* get a column */
2150 col = A [rp++] ;
2151 ASSERT (COL_IS_ALIVE (Col, col) && col != pivot_col) ;
2152 hash = 0 ;
2153 cur_score = 0 ;
2154 cp = Col [col].start ;
2155 /* compact the column */
2156 new_cp = cp ;
2157 cp_end = cp + Col [col].length ;
2159 DEBUG4 ("Adding set diffs for Col: %d.\n", col) ;
2161 while (cp < cp_end)
2163 /* get a row */
2164 row = A [cp++] ;
2165 ASSERT(row >= 0 && row < n_row) ;
2166 row_mark = Row [row].mark ;
2167 /* skip if dead */
2168 if (ROW_IS_MARKED_DEAD (row_mark))
2170 DEBUG4 (" Row %d, dead\n", row) ;
2171 continue ;
2173 DEBUG4 (" Row %d, set diff %d\n", row, row_mark-tag_mark);
2174 ASSERT (row_mark >= tag_mark) ;
2175 /* compact the column */
2176 A [new_cp++] = row ;
2177 /* compute hash function */
2178 hash += row ;
2179 /* add set difference */
2180 cur_score += row_mark - tag_mark ;
2181 /* integer overflow... */
2182 cur_score = MIN (cur_score, n_col) ;
2185 /* recompute the column's length */
2186 Col [col].length = (int) (new_cp - Col [col].start) ;
2188 /* === Further mass elimination ================================= */
2190 if (Col [col].length == 0)
2192 DEBUG4 ("further mass elimination. Col: %d\n", col) ;
2193 /* nothing left but the pivot row in this column */
2194 KILL_PRINCIPAL_COL (Col, col) ;
2195 pivot_row_degree -= Col [col].thickness ;
2196 ASSERT (pivot_row_degree >= 0) ;
2197 /* order it */
2198 Col [col].order = k ;
2199 /* increment order count by column thickness */
2200 k += Col [col].thickness ;
2202 else
2204 /* === Prepare for supercolumn detection ==================== */
2206 DEBUG4 ("Preparing supercol detection for Col: %d.\n", col) ;
2208 /* save score so far */
2209 Col [col].score = cur_score ;
2211 /* add column to hash table, for supercolumn detection */
2212 hash %= n_col + 1 ;
2214 DEBUG4 (" Hash = %d, n_col = %d.\n", hash, n_col) ;
2215 ASSERT (((int) hash) <= n_col) ;
2217 head_column = head [hash] ;
2218 if (head_column > EMPTY)
2220 /* degree list "hash" is non-empty, use prev (shared3) of */
2221 /* first column in degree list as head of hash bucket */
2222 first_col = Col [head_column].headhash ;
2223 Col [head_column].headhash = col ;
2225 else
2227 /* degree list "hash" is empty, use head as hash bucket */
2228 first_col = - (head_column + 2) ;
2229 head [hash] = - (col + 2) ;
2231 Col [col].hash_next = first_col ;
2233 /* save hash function in Col [col].shared3.hash */
2234 Col [col].hash = (int) hash ;
2235 ASSERT (COL_IS_ALIVE (Col, col)) ;
2239 /* The approximate external column degree is now computed. */
2241 /* === Supercolumn detection ======================================== */
2243 DEBUG3 ("** Supercolumn detection phase. **\n") ;
2245 if (!NDEBUG)
2247 detect_super_cols (n_col, Row,
2248 Col, A, head, pivot_row_start, pivot_row_length) ;
2249 } else {
2250 detect_super_cols (Col, A, head, pivot_row_start, pivot_row_length) ;
2253 /* === Kill the pivotal column ====================================== */
2255 KILL_PRINCIPAL_COL (Col, pivot_col) ;
2257 /* === Clear mark =================================================== */
2259 tag_mark = clear_mark (tag_mark+max_deg+1, max_mark, n_row, Row) ;
2261 if (!NDEBUG)
2263 DEBUG3 ("check3\n") ;
2264 debug_mark (n_row, Row, tag_mark, max_mark) ;
2265 } /* NDEBUG */
2267 /* === Finalize the new pivot row, and column scores ================ */
2269 DEBUG3 ("** Finalize scores phase. **\n") ;
2271 /* for each column in pivot row */
2272 rp = pivot_row_start ;
2273 /* compact the pivot row */
2274 new_rp = rp ;
2275 rp_end = rp + pivot_row_length ;
2276 while (rp < rp_end)
2278 col = A [rp++] ;
2279 /* skip dead columns */
2280 if (COL_IS_DEAD (Col, col))
2282 continue ;
2284 A [new_rp++] = col ;
2285 /* add new pivot row to column */
2286 A [Col [col].start + (Col [col].length++)] = pivot_row ;
2288 /* retrieve score so far and add on pivot row's degree. */
2289 /* (we wait until here for this in case the pivot */
2290 /* row's degree was reduced due to mass elimination). */
2291 cur_score = Col [col].score + pivot_row_degree ;
2293 /* calculate the max possible score as the number of */
2294 /* external columns minus the 'k' value minus the */
2295 /* columns thickness */
2296 max_score = n_col - k - Col [col].thickness ;
2298 /* make the score the external degree of the union-of-rows */
2299 cur_score -= Col [col].thickness ;
2301 /* make sure score is less or equal than the max score */
2302 cur_score = MIN (cur_score, max_score) ;
2303 ASSERT (cur_score >= 0) ;
2305 /* store updated score */
2306 Col [col].score = cur_score ;
2308 /* === Place column back in degree list ========================= */
2310 ASSERT (min_score >= 0) ;
2311 ASSERT (min_score <= n_col) ;
2312 ASSERT (cur_score >= 0) ;
2313 ASSERT (cur_score <= n_col) ;
2314 ASSERT (head [cur_score] >= EMPTY) ;
2315 next_col = head [cur_score] ;
2316 Col [col].degree_next = next_col ;
2317 Col [col].prev = EMPTY ;
2318 if (next_col != EMPTY)
2320 Col [next_col].prev = col ;
2322 head [cur_score] = col ;
2324 /* see if this score is less than current min */
2325 min_score = MIN (min_score, cur_score) ;
2329 if (!NDEBUG)
2331 debug_deg_lists (n_row, n_col, Row, Col, head,
2332 min_score, n_col2-k, max_deg) ;
2333 } /* NDEBUG */
2335 /* === Resurrect the new pivot row ================================== */
2337 if (pivot_row_degree > 0)
2339 /* update pivot row length to reflect any cols that were killed */
2340 /* during super-col detection and mass elimination */
2341 Row [pivot_row].start = pivot_row_start ;
2342 Row [pivot_row].length = (int) (new_rp - pivot_row_start) ;
2343 ASSERT (Row [pivot_row].length > 0) ;
2344 Row [pivot_row].degree = pivot_row_degree ;
2345 Row [pivot_row].mark = 0 ;
2346 /* pivot row is no longer dead */
2348 DEBUG1 ("Resurrect Pivot_row %d deg: %d\n",
2349 pivot_row, pivot_row_degree) ;
2353 /* === All principal columns have now been ordered ================== */
2355 return (ngarbage) ;
2359 /* ========================================================================== */
2360 /* === order_children ======================================================= */
2361 /* ========================================================================== */
2364 * The find_ordering routine has ordered all of the principal columns (the
2365 * representatives of the supercolumns). The non-principal columns have not
2366 * yet been ordered. This routine orders those columns by walking up the
2367 * parent tree (a column is a child of the column which absorbed it). The
2368 * final permutation vector is then placed in p [0 ... n_col-1], with p [0]
2369 * being the first column, and p [n_col-1] being the last. It doesn't look
2370 * like it at first glance, but be assured that this routine takes time linear
2371 * in the number of columns. Although not immediately obvious, the time
2372 * taken by this routine is O (n_col), that is, linear in the number of
2373 * columns. Not user-callable.
2375 * @param n_col number of columns of A
2376 * @param Col of size n_col+1
2377 * @param p p [0 ... n_col-1] is the column permutation
2379 private static void order_children (int n_col, Colamd_Col[] Col, int[] p)
2381 /* === Local variables ============================================== */
2383 int i ; /* loop counter for all columns */
2384 int c ; /* column index */
2385 int parent ; /* index of column's parent */
2386 int order ; /* column's order */
2388 /* === Order each non-principal column ============================== */
2390 for (i = 0 ; i < n_col ; i++)
2392 /* find an un-ordered non-principal column */
2393 ASSERT (COL_IS_DEAD (Col, i)) ;
2394 if (!COL_IS_DEAD_PRINCIPAL (Col, i) && Col [i].order == EMPTY)
2396 parent = i ;
2397 /* once found, find its principal parent */
2400 parent = Col [parent].parent ;
2401 } while (!COL_IS_DEAD_PRINCIPAL (Col, parent)) ;
2403 /* now, order all un-ordered non-principal columns along path */
2404 /* to this parent. collapse tree at the same time */
2405 c = i ;
2406 /* get order of parent */
2407 order = Col [parent].order ;
2411 ASSERT (Col [c].order == EMPTY) ;
2413 /* order this column */
2414 Col [c].order = order++ ;
2415 /* collaps tree */
2416 Col [c].parent = parent ;
2418 /* get immediate parent of this column */
2419 c = Col [c].parent ;
2421 /* continue until we hit an ordered column. There are */
2422 /* guarranteed not to be anymore unordered columns */
2423 /* above an ordered column */
2424 } while (Col [c].order == EMPTY) ;
2426 /* re-order the super_col parent to largest order for this group */
2427 Col [parent].order = order ;
2431 /* === Generate the permutation ===================================== */
2433 for (c = 0 ; c < n_col ; c++)
2435 p [Col [c].order] = c ;
2440 /* ========================================================================== */
2441 /* === detect_super_cols ==================================================== */
2442 /* ========================================================================== */
2445 * Detects supercolumns by finding matches between columns in the hash buckets.
2446 * Check amongst columns in the set A [row_start ... row_start + row_length-1].
2447 * The columns under consideration are currently *not* in the degree lists,
2448 * and have already been placed in the hash buckets.
2450 * The hash bucket for columns whose hash function is equal to h is stored
2451 * as follows:
2453 * if head [h] is >= 0, then head [h] contains a degree list, so:
2455 * head [h] is the first column in degree bucket h.
2456 * Col [head [h]].headhash gives the first column in hash bucket h.
2458 * otherwise, the degree list is empty, and:
2460 * -(head [h] + 2) is the first column in hash bucket h.
2462 * For a column c in a hash bucket, Col [c].shared3.prev is NOT a "previous
2463 * column" pointer. Col [c].shared3.hash is used instead as the hash number
2464 * for that column. The value of Col [c].shared4.hash_next is the next column
2465 * in the same hash bucket.
2467 * Assuming no, or "few" hash collisions, the time taken by this routine is
2468 * linear in the sum of the sizes (lengths) of each column whose score has
2469 * just been computed in the approximate degree computation.
2470 * Not user-callable.
2472 * @param Col of size n_col+1
2473 * @param A row indices of A
2474 * @param head head of degree lists and hash buckets
2475 * @param row_start pointer to set of columns to check
2476 * @param row_length number of columns to check
2478 private static void detect_super_cols (Colamd_Col[] Col, int[] A,
2479 int[] head, int row_start, int row_length)
2481 detect_super_cols (0, null, Col, A, head, row_start, row_length) ;
2485 * debug version
2487 * @param n_col number of columns of A
2488 * @param Row of size n_row+1
2489 * @param Col of size n_col+1
2490 * @param A row indices of A
2491 * @param head head of degree lists and hash buckets
2492 * @param row_start pointer to set of columns to check
2493 * @param row_length number of columns to check
2495 private static void detect_super_cols (int n_col, Colamd_Row[] Row,
2496 Colamd_Col[] Col, int[] A, int[] head, int row_start, int row_length)
2498 /* === Local variables ============================================== */
2500 int hash ; /* hash value for a column */
2501 int rp ; /* pointer to a row */
2502 int c ; /* a column index */
2503 int super_c ; /* column index of the column to absorb into */
2504 int cp1 ; /* column pointer for column super_c */
2505 int cp2 ; /* column pointer for column c */
2506 int length ; /* length of column super_c */
2507 int prev_c ; /* column preceding c in hash bucket */
2508 int i ; /* loop counter */
2509 int rp_end ; /* pointer to the end of the row */
2510 int col ; /* a column index in the row to check */
2511 int head_column ; /* first column in hash bucket or degree list */
2512 int first_col ; /* first column in hash bucket */
2514 /* === Consider each column in the row ============================== */
2516 rp = row_start ;
2517 rp_end = rp + row_length ;
2518 while (rp < rp_end)
2520 col = A [rp++] ;
2521 if (COL_IS_DEAD (Col, col))
2523 continue ;
2526 /* get hash number for this column */
2527 hash = Col [col].hash ;
2528 ASSERT (hash <= n_col) ;
2530 /* === Get the first column in this hash bucket ===================== */
2532 head_column = head [hash] ;
2533 if (head_column > EMPTY)
2535 first_col = Col [head_column].headhash ;
2537 else
2539 first_col = - (head_column + 2) ;
2542 /* === Consider each column in the hash bucket ====================== */
2544 for (super_c = first_col ; super_c != EMPTY ;
2545 super_c = Col [super_c].hash_next)
2547 ASSERT (COL_IS_ALIVE (Col, super_c)) ;
2548 ASSERT (Col [super_c].hash == hash) ;
2549 length = Col [super_c].length ;
2551 /* prev_c is the column preceding column c in the hash bucket */
2552 prev_c = super_c ;
2554 /* === Compare super_c with all columns after it ================ */
2556 for (c = Col [super_c].hash_next ;
2557 c != EMPTY ; c = Col [c].hash_next)
2559 ASSERT (c != super_c) ;
2560 ASSERT (COL_IS_ALIVE (Col, c)) ;
2561 ASSERT (Col [c].hash == hash) ;
2563 /* not identical if lengths or scores are different */
2564 if (Col [c].length != length ||
2565 Col [c].score != Col [super_c].score)
2567 prev_c = c ;
2568 continue ;
2571 /* compare the two columns */
2572 cp1 = Col [super_c].start ;
2573 cp2 = Col [c].start ;
2575 for (i = 0 ; i < length ; i++)
2577 /* the columns are "clean" (no dead rows) */
2578 ASSERT (ROW_IS_ALIVE (Row, A [cp1])) ;
2579 ASSERT (ROW_IS_ALIVE (Row, A [cp2])) ;
2580 /* row indices will same order for both supercols, */
2581 /* no gather scatter nessasary */
2582 if (A [cp1++] != A [cp2++])
2584 break ;
2588 /* the two columns are different if the for-loop "broke" */
2589 if (i != length)
2591 prev_c = c ;
2592 continue ;
2595 /* === Got it! two columns are identical =================== */
2597 ASSERT (Col [c].score == Col [super_c].score) ;
2599 Col [super_c].thickness += Col [c].thickness ;
2600 Col [c].parent = super_c ;
2601 KILL_NON_PRINCIPAL_COL (Col, c) ;
2602 /* order c later, in order_children() */
2603 Col [c].order = EMPTY ;
2604 /* remove c from hash bucket */
2605 Col [prev_c].hash_next = Col [c].hash_next ;
2609 /* === Empty this hash bucket ======================================= */
2611 if (head_column > EMPTY)
2613 /* corresponding degree list "hash" is not empty */
2614 Col [head_column].headhash = EMPTY ;
2616 else
2618 /* corresponding degree list "hash" is empty */
2619 head [hash] = EMPTY ;
2625 /* ========================================================================== */
2626 /* === garbage_collection =================================================== */
2627 /* ========================================================================== */
2630 * Defragments and compacts columns and rows in the workspace A. Used when
2631 * all avaliable memory has been used while performing row merging. Returns
2632 * the index of the first free position in A, after garbage collection. The
2633 * time taken by this routine is linear in the size of the array A, which is
2634 * itself linear in the number of nonzeros in the input matrix.
2635 * Not user-callable.
2637 * @param n_row number of rows
2638 * @param n_col number of columns
2639 * @param Row row info
2640 * @param Col column info
2641 * @param A A [0 ... Alen-1] holds the matrix
2642 * @param pfree
2643 * @return the new value of pfree
2645 private static int garbage_collection (int n_row, int n_col,
2646 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int pfree)
2648 /* === Local variables ============================================== */
2650 int psrc ; /* source pointer */
2651 int pdest ; /* destination pointer */
2652 int j ; /* counter */
2653 int r ; /* a row index */
2654 int c ; /* a column index */
2655 int length ; /* length of a row or column */
2657 int debug_rows = 0 ;
2658 if (!NDEBUG)
2660 DEBUG2 ("Defrag..\n") ;
2661 for (psrc = 0 ; psrc < pfree ; psrc++) ASSERT (A [psrc] >= 0) ;
2662 debug_rows = 0 ;
2665 /* === Defragment the columns ======================================= */
2667 pdest = 0 ;
2668 for (c = 0 ; c < n_col ; c++)
2670 if (COL_IS_ALIVE (Col, c))
2672 psrc = Col [c].start ;
2674 /* move and compact the column */
2675 ASSERT (pdest <= psrc) ;
2676 Col [c].start = (int) (pdest - 0) ;
2677 length = Col [c].length ;
2678 for (j = 0 ; j < length ; j++)
2680 r = A [psrc++] ;
2681 if (ROW_IS_ALIVE (Row, r))
2683 A [pdest] = r ;
2686 Col [c].length = (int) (pdest - Col [c].start) ;
2690 /* === Prepare to defragment the rows =============================== */
2692 for (r = 0 ; r < n_row ; r++)
2694 if (ROW_IS_DEAD (Row, r) || (Row [r].length == 0))
2696 /* This row is already dead, or is of zero length. Cannot compact
2697 * a row of zero length, so kill it. NOTE: in the current version,
2698 * there are no zero-length live rows. Kill the row (for the first
2699 * time, or again) just to be safe. */
2700 KILL_ROW (Row, r) ;
2702 else
2704 /* save first column index in Row [r].shared2.first_column */
2705 psrc = Row [r].start ;
2706 Row [r].first_column = A [psrc] ;
2707 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2708 /* flag the start of the row with the one's complement of row */
2709 A [psrc] = ONES_COMPLEMENT (r) ;
2710 if (!NDEBUG)
2712 debug_rows++ ;
2713 } /* NDEBUG */
2717 /* === Defragment the rows ========================================== */
2719 psrc = pdest ;
2720 while (psrc < pfree)
2722 /* find a negative number ... the start of a row */
2723 if (A [psrc++] < 0)
2725 psrc-- ;
2726 /* get the row index */
2727 r = ONES_COMPLEMENT (A [psrc]) ;
2728 ASSERT (r >= 0 && r < n_row) ;
2729 /* restore first column index */
2730 A [psrc] = Row [r].first_column ;
2731 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2732 ASSERT (Row [r].length > 0) ;
2733 /* move and compact the row */
2734 ASSERT (pdest <= psrc) ;
2735 Row [r].start = (int) (pdest - 0) ;
2736 length = Row [r].length ;
2737 for (j = 0 ; j < length ; j++)
2739 c = A [psrc++] ;
2740 if (COL_IS_ALIVE (Col, c))
2742 A [pdest++] = c ;
2745 Row [r].length = (int) (pdest - Row [r].start) ;
2746 ASSERT (Row [r].length > 0) ;
2747 if (!NDEBUG)
2749 debug_rows-- ;
2750 } /* NDEBUG */
2753 /* ensure we found all the rows */
2754 ASSERT (debug_rows == 0) ;
2756 /* === Return the new value of pfree ================================ */
2758 return ((int) (pdest - 0)) ;
2762 /* ========================================================================== */
2763 /* === clear_mark =========================================================== */
2764 /* ========================================================================== */
2767 * Clears the Row [].shared2.mark array, and returns the new tag_mark.
2769 * @param tag_mark new value of tag_mark
2770 * @param max_mark max allowed value of tag_mark
2771 * @param n_row number of rows in A
2772 * @param Row Row [0 ... n_row-1].shared2.mark is set to zero
2773 * @return the new value for tag_mark
2775 private static int clear_mark (int tag_mark, int max_mark, int n_row,
2776 Colamd_Row[] Row)
2778 int r ;
2780 if (tag_mark <= 0 || tag_mark >= max_mark)
2782 for (r = 0 ; r < n_row ; r++)
2784 if (ROW_IS_ALIVE (Row, r))
2786 Row [r].mark = 0 ;
2789 tag_mark = 1 ;
2792 return (tag_mark) ;
2796 /* ========================================================================== */
2797 /* === print_report ========================================================= */
2798 /* ========================================================================== */
2802 * @param method
2803 * @param stats
2805 private static void print_report (String method, int[] stats)
2807 int i1, i2, i3 ;
2809 PRINTF ("\n%s version %d.%d, %s: ", method,
2810 COLAMD_MAIN_VERSION, COLAMD_SUB_VERSION, COLAMD_DATE) ;
2812 if (stats == null)
2814 PRINTF ("No statistics available.\n") ;
2815 return ;
2818 i1 = stats [COLAMD_INFO1] ;
2819 i2 = stats [COLAMD_INFO2] ;
2820 i3 = stats [COLAMD_INFO3] ;
2822 if (stats [COLAMD_STATUS] >= 0)
2824 PRINTF ("OK. ") ;
2826 else
2828 PRINTF ("ERROR. ") ;
2831 switch (stats [COLAMD_STATUS])
2834 case COLAMD_OK_BUT_JUMBLED:
2836 PRINTF("Matrix has unsorted or duplicate row indices.\n") ;
2838 PRINTF("%s: number of duplicate or out-of-order row indices: %d\n",
2839 method, i3) ;
2841 PRINTF("%s: last seen duplicate or out-of-order row index: %d\n",
2842 method, INDEX (i2)) ;
2844 PRINTF("%s: last seen in column: %d",
2845 method, INDEX (i1)) ;
2847 /* no break - fall through to next case instead */
2849 case COLAMD_OK:
2851 PRINTF("\n") ;
2853 PRINTF("%s: number of dense or empty rows ignored: %d\n",
2854 method, stats [COLAMD_DENSE_ROW]) ;
2856 PRINTF("%s: number of dense or empty columns ignored: %d\n",
2857 method, stats [COLAMD_DENSE_COL]) ;
2859 PRINTF("%s: number of garbage collections performed: %d\n",
2860 method, stats [COLAMD_DEFRAG_COUNT]) ;
2861 break ;
2863 case COLAMD_ERROR_A_not_present:
2865 PRINTF("Array A (row indices of matrix) not present.\n") ;
2866 break ;
2868 case COLAMD_ERROR_p_not_present:
2870 PRINTF("Array p (column pointers for matrix) not present.\n") ;
2871 break ;
2873 case COLAMD_ERROR_nrow_negative:
2875 PRINTF("Invalid number of rows (%d).\n", i1) ;
2876 break ;
2878 case COLAMD_ERROR_ncol_negative:
2880 PRINTF("Invalid number of columns (%d).\n", i1) ;
2881 break ;
2883 case COLAMD_ERROR_nnz_negative:
2885 PRINTF("Invalid number of nonzero entries (%d).\n", i1) ;
2886 break ;
2888 case COLAMD_ERROR_p0_nonzero:
2890 PRINTF("Invalid column pointer, p [0] = %d, must be zero.\n", i1);
2891 break ;
2893 case COLAMD_ERROR_A_too_small:
2895 PRINTF("Array A too small.\n") ;
2896 PRINTF(" Need Alen >= %d, but given only Alen = %d.\n",
2897 i1, i2) ;
2898 break ;
2900 case COLAMD_ERROR_col_length_negative:
2902 PRINTF
2903 ("Column %d has a negative number of nonzero entries (%d).\n",
2904 INDEX (i1), i2) ;
2905 break ;
2907 case COLAMD_ERROR_row_index_out_of_bounds:
2909 PRINTF
2910 ("Row index (row %d) out of bounds (%d to %d) in column %d.\n",
2911 INDEX (i2), INDEX (0), INDEX (i3-1), INDEX (i1)) ;
2912 break ;
2914 case COLAMD_ERROR_out_of_memory:
2916 PRINTF("Out of memory.\n") ;
2917 break ;
2919 /* v2.4: internal-error case deleted */
2924 /* ========================================================================== */
2925 /* === colamd debugging routines ============================================ */
2926 /* ========================================================================== */
2929 * At this point, all empty rows and columns are dead. All live columns
2930 * are "clean" (containing no dead rows) and simplicial (no supercolumns
2931 * yet). Rows may contain dead columns, but all live rows contain at
2932 * least one live column.
2934 * @param n_row
2935 * @param n_col
2936 * @param Row
2937 * @param Col
2938 * @param A
2939 * @param n_col2
2941 private static void debug_structures (int n_row, int n_col,
2942 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int n_col2)
2944 if (!NDEBUG)
2946 /* === Local variables ============================================== */
2948 int i ;
2949 int c ;
2950 int cp ;
2951 int cp_end ;
2952 int len ;
2953 int score ;
2954 int r ;
2955 int rp ;
2956 int rp_end ;
2957 int deg ;
2959 /* === Check A, Row, and Col ======================================== */
2961 for (c = 0 ; c < n_col ; c++)
2963 if (COL_IS_ALIVE (Col, c))
2965 len = Col [c].length ;
2966 score = Col [c].score ;
2967 DEBUG4 ("initial live col %5d %5d %5d\n", c, len, score) ;
2968 ASSERT (len > 0) ;
2969 ASSERT (score >= 0) ;
2970 ASSERT (Col [c].thickness == 1) ;
2971 cp = Col [c].start ;
2972 cp_end = cp + len ;
2973 while (cp < cp_end)
2975 r = A [cp++] ;
2976 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2979 else
2981 i = Col [c].order ;
2982 ASSERT (i >= n_col2 && i < n_col) ;
2986 for (r = 0 ; r < n_row ; r++)
2988 if (ROW_IS_ALIVE (Row, r))
2990 i = 0 ;
2991 len = Row [r].length ;
2992 deg = Row [r].degree ;
2993 ASSERT (len > 0) ;
2994 ASSERT (deg > 0) ;
2995 rp = Row [r].start ;
2996 rp_end = rp + len ;
2997 while (rp < rp_end)
2999 c = A [rp++] ;
3000 if (COL_IS_ALIVE (Col, c))
3002 i++ ;
3005 ASSERT (i > 0) ;
3012 /* ========================================================================== */
3013 /* === debug_deg_lists ====================================================== */
3014 /* ========================================================================== */
3017 * Prints the contents of the degree lists. Counts the number of columns
3018 * in the degree list and compares it to the total it should have. Also
3019 * checks the row degrees.
3021 * @param n_row
3022 * @param n_col
3023 * @param Row
3024 * @param Col
3025 * @param head
3026 * @param min_score
3027 * @param should
3028 * @param max_deg
3030 private static void debug_deg_lists (int n_row, int n_col,
3031 Colamd_Row[] Row, Colamd_Col[] Col, int[] head, int min_score,
3032 int should, int max_deg)
3034 if (!NDEBUG)
3036 /* === Local variables ============================================== */
3038 int deg ;
3039 int col ;
3040 int have ;
3041 int row ;
3043 /* === Check the degree lists ======================================= */
3045 if (n_col > 10000 && colamd_debug <= 0)
3047 return ;
3049 have = 0 ;
3050 DEBUG4 ("Degree lists: %d\n", min_score) ;
3051 for (deg = 0 ; deg <= n_col ; deg++)
3053 col = head [deg] ;
3054 if (col == EMPTY)
3056 continue ;
3058 DEBUG4 ("%d:", deg) ;
3059 while (col != EMPTY)
3061 DEBUG4 (" %d", col) ;
3062 have += Col [col].thickness ;
3063 ASSERT (COL_IS_ALIVE (Col, col)) ;
3064 col = Col [col].degree_next ;
3066 DEBUG4 ("\n") ;
3068 DEBUG4 ("should %d have %d\n", should, have) ;
3069 ASSERT (should == have) ;
3071 /* === Check the row degrees ======================================== */
3073 if (n_row > 10000 && colamd_debug <= 0)
3075 return ;
3077 for (row = 0 ; row < n_row ; row++)
3079 if (ROW_IS_ALIVE (Row, row))
3081 ASSERT (Row [row].degree <= max_deg) ;
3088 /* ========================================================================== */
3089 /* === debug_mark =========================================================== */
3090 /* ========================================================================== */
3093 * Ensures that the tag_mark is less that the maximum and also ensures that
3094 * each entry in the mark array is less than the tag mark.
3096 * @param n_row
3097 * @param Row
3098 * @param tag_mark
3099 * @param max_mark
3101 private static void debug_mark (int n_row, Colamd_Row[] Row, int tag_mark,
3102 int max_mark)
3104 if (!NDEBUG)
3106 /* === Local variables ============================================== */
3108 int r ;
3110 /* === Check the Row marks ========================================== */
3112 ASSERT (tag_mark > 0 && tag_mark <= max_mark) ;
3113 if (n_row > 10000 && colamd_debug <= 0)
3115 return ;
3117 for (r = 0 ; r < n_row ; r++)
3119 ASSERT (Row [r].mark < tag_mark) ;
3125 /* ========================================================================== */
3126 /* === debug_matrix ========================================================= */
3127 /* ========================================================================== */
3130 * Prints out the contents of the columns and the rows.
3132 * @param n_row
3133 * @param n_col
3134 * @param Row
3135 * @param Col
3136 * @param A
3138 private static void debug_matrix (int n_row, int n_col,
3139 Colamd_Row[] Row, Colamd_Col[] Col, int[] A)
3141 if (!NDEBUG)
3143 /* === Local variables ============================================== */
3145 int r ;
3146 int c ;
3147 int rp ;
3148 int rp_end ;
3149 int cp ;
3150 int cp_end ;
3152 /* === Dump the rows and columns of the matrix ====================== */
3154 if (colamd_debug < 3)
3156 return ;
3158 DEBUG3 ("DUMP MATRIX:\n") ;
3159 for (r = 0 ; r < n_row ; r++)
3161 DEBUG3 ("Row %d alive? %d\n", r, ROW_IS_ALIVE (Row, r) ? 1 : 0) ;
3162 if (ROW_IS_DEAD (Row, r))
3164 continue ;
3166 DEBUG3 ("start %d length %d degree %d\n",
3167 Row [r].start, Row [r].length, Row [r].degree) ;
3168 rp = Row [r].start ;
3169 rp_end = rp + Row [r].length ;
3170 while (rp < rp_end)
3172 c = A [rp++] ;
3173 DEBUG4 (" %d col %d\n", COL_IS_ALIVE (Col, c) ? 1 : 0, c) ;
3177 for (c = 0 ; c < n_col ; c++)
3179 DEBUG3 ("Col %d alive? %d\n", c, COL_IS_ALIVE (Col, c) ? 1 : 0) ;
3180 if (COL_IS_DEAD (Col, c))
3182 continue ;
3184 DEBUG3 ("start %d length %d shared1 %d shared2 %d\n",
3185 Col [c].start, Col [c].length,
3186 Col [c].thickness, Col [c].score) ;
3187 cp = Col [c].start ;
3188 cp_end = cp + Col [c].length ;
3189 while (cp < cp_end)
3191 r = A [cp++] ;
3192 DEBUG4 (" %d row %d\n", ROW_IS_ALIVE (Row, r) ? 1 : 0, r) ;
3198 @SuppressWarnings("unused")
3199 private static void colamd_get_debug (String method)
3201 if (!NDEBUG)
3203 File f ;
3204 colamd_debug = 0 ; /* no debug printing */
3205 f = new File("debug") ;
3206 if (!f.exists())
3208 colamd_debug = 0 ;
3210 else
3212 try {
3213 FileReader fr ;
3214 fr = new FileReader(f) ;
3215 BufferedReader br ;
3216 br = new BufferedReader(fr) ;
3217 colamd_debug = Integer.valueOf( br.readLine() ) ;
3218 br.close() ;
3219 fr.close() ;
3220 } catch (IOException e) {
3221 PRINTF ("%s: AMD_debug_init, " +
3222 "error reading debug.amd file", method) ;
3225 DEBUG0 ("%s: debug version, D = %d (THIS WILL BE SLOW!)\n",
3226 method, colamd_debug) ;