Squashing warnings.
[COLAMDJ.git] / src / main / java / edu / ufl / cise / colamd / tdouble / Dcolamd.java
blob980d45a3f37e46bd553edfe403ef5ff3f8406e48
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 ;
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 col ; /* a column index */
1307 int row ; /* a row index */
1308 int cp ; /* a column pointer */
1309 int cp_end ; /* a pointer to the end of a column */
1310 int rp ; /* a row pointer */
1311 int rp_end ; /* a pointer to the end of a row */
1312 int last_row ; /* previous row */
1314 /* === Initialize columns, and check column pointers ================ */
1316 for (col = 0 ; col < n_col ; col++)
1318 Col [col].start = p [col] ;
1319 Col [col].length = p [col+1] - p [col] ;
1321 if (Col [col].length < 0)
1323 /* column pointers must be non-decreasing */
1324 stats [COLAMD_STATUS] = COLAMD_ERROR_col_length_negative ;
1325 stats [COLAMD_INFO1] = col ;
1326 stats [COLAMD_INFO2] = Col [col].length ;
1327 DEBUG0 ("colamd: col %d length %d < 0\n", col, Col [col].length) ;
1328 return (FALSE) ;
1331 Col [col].thickness(1) ;
1332 Col [col].score(0) ;
1333 Col [col].prev(EMPTY) ;
1334 Col [col].degree_next(EMPTY) ;
1337 /* p [0..n_col] no longer needed, used as "head" in subsequent routines */
1339 /* === Scan columns, compute row degrees, and check row indices ===== */
1341 stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
1343 for (row = 0 ; row < n_row ; row++)
1345 Row [row].length = 0 ;
1346 Row [row].mark(-1) ;
1349 for (col = 0 ; col < n_col ; col++)
1351 last_row = -1 ;
1353 cp = p [col] ;
1354 cp_end = p [col+1] ;
1356 while (cp < cp_end)
1358 row = A [cp++] ;
1360 /* make sure row indices within range */
1361 if (row < 0 || row >= n_row)
1363 stats [COLAMD_STATUS] = COLAMD_ERROR_row_index_out_of_bounds ;
1364 stats [COLAMD_INFO1] = col ;
1365 stats [COLAMD_INFO2] = row ;
1366 stats [COLAMD_INFO3] = n_row ;
1367 DEBUG0 ("colamd: row %d col %d out of bounds\n", row, col) ;
1368 return (FALSE) ;
1371 if (row <= last_row || Row [row].mark() == col)
1373 /* row index are unsorted or repeated (or both), thus col */
1374 /* is jumbled. This is a notice, not an error condition. */
1375 stats [COLAMD_STATUS] = COLAMD_OK_BUT_JUMBLED ;
1376 stats [COLAMD_INFO1] = col ;
1377 stats [COLAMD_INFO2] = row ;
1378 (stats [COLAMD_INFO3]) ++ ;
1379 DEBUG1 ("colamd: row %d col %d unsorted/duplicate\n",row,col);
1382 if (Row [row].mark() != col)
1384 Row [row].length++ ;
1386 else
1388 /* this is a repeated entry in the column, */
1389 /* it will be removed */
1390 Col [col].length-- ;
1393 /* mark the row as having been seen in this column */
1394 Row [row].mark(col) ;
1396 last_row = row ;
1400 /* === Compute row pointers ========================================= */
1402 /* row form of the matrix starts directly after the column */
1403 /* form of matrix in A */
1404 Row [0].start = p [n_col] ;
1405 Row [0].p(Row [0].start) ;
1406 Row [0].mark(-1) ;
1407 for (row = 1 ; row < n_row ; row++)
1409 Row [row].start = Row [row-1].start + Row [row-1].length ;
1410 Row [row].p(Row [row].start) ;
1411 Row [row].mark(-1) ;
1414 /* === Create row form ============================================== */
1416 if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1418 /* if cols jumbled, watch for repeated row indices */
1419 for (col = 0 ; col < n_col ; col++)
1421 cp = p [col] ;
1422 cp_end = p [col+1] ;
1424 while (cp < cp_end)
1426 row = A [cp++] ;
1428 if (Row [row].mark() != col)
1430 A [Row [row].p()] = col ;
1431 Row [row].p( Row [row].p() + 1 ) ;
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 ;
1447 Row [A [cp]].p( Row [A [cp]].p() + 1 ) ;
1448 cp++ ;
1453 /* === Clear the row marks and set row degrees ====================== */
1455 for (row = 0 ; row < n_row ; row++)
1457 Row [row].mark(0) ;
1458 Row [row].degree(Row [row].length) ;
1461 /* === See if we need to re-create columns ========================== */
1463 if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1465 DEBUG0 ("colamd: reconstructing column form, matrix jumbled\n") ;
1467 if (!NDEBUG)
1469 /* make sure column lengths are correct */
1470 for (col = 0 ; col < n_col ; col++)
1472 p [col] = Col [col].length ;
1474 for (row = 0 ; row < n_row ; row++)
1476 rp = Row [row].start ;
1477 rp_end = rp + Row [row].length ;
1478 while (rp < rp_end)
1480 p [A [rp++]]-- ;
1483 for (col = 0 ; col < n_col ; col++)
1485 ASSERT (p [col] == 0) ;
1487 /* now p is all zero (different than when debugging is turned off) */
1488 } /* NDEBUG */
1490 /* === Compute col pointers ========================================= */
1492 /* col form of the matrix starts at A [0]. */
1493 /* Note, we may have a gap between the col form and the row */
1494 /* form if there were duplicate entries, if so, it will be */
1495 /* removed upon the first garbage collection */
1496 Col [0].start = 0 ;
1497 p [0] = Col [0].start ;
1498 for (col = 1 ; col < n_col ; col++)
1500 /* note that the lengths here are for pruned columns, i.e. */
1501 /* no duplicate row indices will exist for these columns */
1502 Col [col].start = Col [col-1].start + Col [col-1].length ;
1503 p [col] = Col [col].start ;
1506 /* === Re-create col form =========================================== */
1508 for (row = 0 ; row < n_row ; row++)
1510 rp = Row [row].start ;
1511 rp_end = rp + Row [row].length ;
1512 while (rp < rp_end)
1514 A [(p [A [rp++]])++] = row ;
1519 /* === Done. Matrix is not (or no longer) jumbled ================== */
1521 return (TRUE) ;
1525 /* ========================================================================== */
1526 /* === init_scoring ========================================================= */
1527 /* ========================================================================== */
1530 * Kills dense or empty columns and rows, calculates an initial score for
1531 * each column, and places all columns in the degree lists.init_rows_cols
1533 * @param n_row number of rows of A
1534 * @param n_col number of columns of A
1535 * @param Row of size n_row+1
1536 * @param Col of size n_col+1
1537 * @param A column form and row form of A
1538 * @param head of size n_col+1
1539 * @param knobs parameters
1540 * @param p_n_row2 size 1, number of non-dense, non-empty rows
1541 * @param p_n_col2 size 1, number of non-dense, non-empty columns
1542 * @param p_max_deg size 1, maximum row degree
1544 private static void init_scoring (int n_row, int n_col, Colamd_Row[] Row,
1545 Colamd_Col[] Col, int[] A, int[] head, double[] knobs,
1546 int[] p_n_row2, int[] p_n_col2, int[] p_max_deg)
1548 /* === Local variables ============================================== */
1550 int c ; /* a column index */
1551 int r, row ; /* a row index */
1552 int cp ; /* a column pointer */
1553 int deg ; /* degree of a row or column */
1554 int cp_end ; /* a pointer to the end of a column */
1555 int new_cp ; /* new column pointer */
1556 int col_length ; /* length of pruned column */
1557 int score ; /* current column score */
1558 int n_col2 ; /* number of non-dense, non-empty columns */
1559 int n_row2 ; /* number of non-dense, non-empty rows */
1560 int dense_row_count ; /* remove rows with more entries than this */
1561 int dense_col_count ; /* remove cols with more entries than this */
1562 int min_score ; /* smallest column score */
1563 int max_deg ; /* maximum row degree */
1564 int next_col ; /* Used to add to degree list.*/
1566 int debug_count = 0 ; /* debug only. */
1568 /* === Extract knobs ================================================ */
1570 /* Note: if knobs contains a NaN, this is undefined: */
1571 if (knobs [COLAMD_DENSE_ROW] < 0)
1573 /* only remove completely dense rows */
1574 dense_row_count = n_col-1 ;
1576 else
1578 dense_row_count = DENSE_DEGREE (knobs [COLAMD_DENSE_ROW], n_col) ;
1580 if (knobs [COLAMD_DENSE_COL] < 0)
1582 /* only remove completely dense columns */
1583 dense_col_count = n_row-1 ;
1585 else
1587 dense_col_count = DENSE_DEGREE (knobs [COLAMD_DENSE_COL],
1588 MIN (n_row, n_col)) ;
1591 DEBUG1 ("colamd: densecount: %d %d\n", dense_row_count, dense_col_count) ;
1592 max_deg = 0 ;
1593 n_col2 = n_col ;
1594 n_row2 = n_row ;
1596 /* === Kill empty columns =========================================== */
1598 /* Put the empty columns at the end in their natural order, so that LU */
1599 /* factorization can proceed as far as possible. */
1600 for (c = n_col-1 ; c >= 0 ; c--)
1602 deg = Col [c].length ;
1603 if (deg == 0)
1605 /* this is a empty column, kill and order it last */
1606 Col [c].order(--n_col2) ;
1607 KILL_PRINCIPAL_COL (Col, c) ;
1610 DEBUG1 ("colamd: null columns killed: %d\n", n_col - n_col2) ;
1612 /* === Kill dense columns =========================================== */
1614 /* Put the dense columns at the end, in their natural order */
1615 for (c = n_col-1 ; c >= 0 ; c--)
1617 /* skip any dead columns */
1618 if (COL_IS_DEAD (Col, c))
1620 continue ;
1622 deg = Col [c].length ;
1623 if (deg > dense_col_count)
1625 /* this is a dense column, kill and order it last */
1626 Col [c].order(--n_col2) ;
1627 /* decrement the row degrees */
1628 cp = Col [c].start ;
1629 cp_end = cp + Col [c].length ;
1630 while (cp < cp_end)
1632 Row [A [cp]].degree( Row [A [cp]].degree() - 1 ) ;
1633 cp++ ;
1635 KILL_PRINCIPAL_COL (Col, c) ;
1638 DEBUG1 ("colamd: Dense and null columns killed: %d\n", n_col - n_col2) ;
1640 /* === Kill dense and empty rows ==================================== */
1642 for (r = 0 ; r < n_row ; r++)
1644 deg = Row [r].degree() ;
1645 ASSERT (deg >= 0 && deg <= n_col) ;
1646 if (deg > dense_row_count || deg == 0)
1648 /* kill a dense or empty row */
1649 KILL_ROW (Row, r) ;
1650 --n_row2 ;
1652 else
1654 /* keep track of max degree of remaining rows */
1655 max_deg = MAX (max_deg, deg) ;
1658 DEBUG1 ("colamd: Dense and null rows killed: %d\n", n_row - n_row2) ;
1660 /* === Compute initial column scores ================================ */
1662 /* At this point the row degrees are accurate. They reflect the number */
1663 /* of "live" (non-dense) columns in each row. No empty rows exist. */
1664 /* Some "live" columns may contain only dead rows, however. These are */
1665 /* pruned in the code below. */
1667 /* now find the initial matlab score for each column */
1668 for (c = n_col-1 ; c >= 0 ; c--)
1670 /* skip dead column */
1671 if (COL_IS_DEAD (Col, c))
1673 continue ;
1675 score = 0 ;
1676 cp = Col [c].start ;
1677 new_cp = cp ;
1678 cp_end = cp + Col [c].length ;
1679 while (cp < cp_end)
1681 /* get a row */
1682 row = A [cp++] ;
1683 /* skip if dead */
1684 if (ROW_IS_DEAD (Row, row))
1686 continue ;
1688 /* compact the column */
1689 A [new_cp++] = row ;
1690 /* add row's external degree */
1691 score += Row [row].degree() - 1 ;
1692 /* guard against integer overflow */
1693 score = MIN (score, n_col) ;
1695 /* determine pruned column length */
1696 col_length = (new_cp - Col [c].start) ;
1697 if (col_length == 0)
1699 /* a newly-made null column (all rows in this col are "dense" */
1700 /* and have already been killed) */
1701 DEBUG2 ("Newly null killed: %d\n", c) ;
1702 Col [c].order(--n_col2) ;
1703 KILL_PRINCIPAL_COL (Col, c) ;
1705 else
1707 /* set column length and set score */
1708 ASSERT (score >= 0) ;
1709 ASSERT (score <= n_col) ;
1710 Col [c].length = col_length ;
1711 Col [c].score(score) ;
1714 DEBUG1 ("colamd: Dense, null, and newly-null columns killed: %d\n",
1715 n_col-n_col2) ;
1717 /* At this point, all empty rows and columns are dead. All live columns */
1718 /* are "clean" (containing no dead rows) and simplicial (no supercolumns */
1719 /* yet). Rows may contain dead columns, but all live rows contain at */
1720 /* least one live column. */
1722 if (!NDEBUG)
1724 debug_structures (n_row, n_col, Row, Col, A, n_col2) ;
1727 /* === Initialize degree lists ========================================== */
1729 if (!NDEBUG)
1731 debug_count = 0 ;
1734 /* clear the hash buckets */
1735 for (c = 0 ; c <= n_col ; c++)
1737 head [c] = EMPTY ;
1739 min_score = n_col ;
1740 /* place in reverse order, so low column indices are at the front */
1741 /* of the lists. This is to encourage natural tie-breaking */
1742 for (c = n_col-1 ; c >= 0 ; c--)
1744 /* only add principal columns to degree lists */
1745 if (COL_IS_ALIVE (Col, c))
1747 DEBUG4 ("place %d score %d minscore %d ncol %d\n",
1748 c, Col [c].score(), min_score, n_col) ;
1750 /* === Add columns score to DList =============================== */
1752 score = Col [c].score() ;
1754 ASSERT (min_score >= 0) ;
1755 ASSERT (min_score <= n_col) ;
1756 ASSERT (score >= 0) ;
1757 ASSERT (score <= n_col) ;
1758 ASSERT (head [score] >= EMPTY) ;
1760 /* now add this column to dList at proper score location */
1761 next_col = head [score] ;
1762 Col [c].prev(EMPTY) ;
1763 Col [c].degree_next(next_col) ;
1765 /* if there already was a column with the same score, set its */
1766 /* previous pointer to this new column */
1767 if (next_col != EMPTY)
1769 Col [next_col].prev(c) ;
1771 head [score] = c ;
1773 /* see if this score is less than current min */
1774 min_score = MIN (min_score, score) ;
1776 if (!NDEBUG)
1778 debug_count++ ;
1784 if (!NDEBUG)
1786 DEBUG1 ("colamd: Live cols %d out of %d, non-princ: %d\n",
1787 debug_count, n_col, n_col-debug_count) ;
1788 ASSERT (debug_count == n_col2) ;
1789 debug_deg_lists (n_row, n_col, Row, Col, head, min_score, n_col2, max_deg) ;
1790 } /* NDEBUG */
1792 /* === Return number of remaining columns, and max row degree ======= */
1794 p_n_col2[0] = n_col2 ;
1795 p_n_row2[0] = n_row2 ;
1796 p_max_deg[0] = max_deg ;
1800 /* ========================================================================== */
1801 /* === find_ordering ======================================================== */
1802 /* ========================================================================== */
1805 * Order the principal columns of the supercolumn form of the matrix
1806 * (no supercolumns on input). Uses a minimum approximate column minimum
1807 * degree ordering method. Not user-callable.
1809 * @param n_row number of rows of A
1810 * @param n_col number of columns of A
1811 * @param Alen size of A, 2*nnz + n_col or larger
1812 * @param Row of size n_row+1
1813 * @param Col of size n_col+1
1814 * @param A column form and row form of A
1815 * @param head of size n_col+1
1816 * @param n_col2 Remaining columns to order
1817 * @param max_deg Maximum row degree
1818 * @param pfree index of first free slot (2*nnz on entry)
1819 * @param aggressive
1820 * @return the number of garbage collections
1822 private static int find_ordering (int n_row, int n_col, int Alen,
1823 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int[] head, int n_col2,
1824 int max_deg, int pfree, int aggressive)
1826 /* === Local variables ============================================== */
1828 int k ; /* current pivot ordering step */
1829 int pivot_col ; /* current pivot column */
1830 int cp ; /* a column pointer */
1831 int rp ; /* a row pointer */
1832 int pivot_row ; /* current pivot row */
1833 int new_cp ; /* modified column pointer */
1834 int new_rp ; /* modified row pointer */
1835 int pivot_row_start ; /* pointer to start of pivot row */
1836 int pivot_row_degree ; /* number of columns in pivot row */
1837 int pivot_row_length ; /* number of supercolumns in pivot row */
1838 int pivot_col_score ; /* score of pivot column */
1839 int needed_memory ; /* free space needed for pivot row */
1840 int cp_end ; /* pointer to the end of a column */
1841 int rp_end ; /* pointer to the end of a row */
1842 int row ; /* a row index */
1843 int col ; /* a column index */
1844 int max_score ; /* maximum possible score */
1845 int cur_score ; /* score of current column */
1846 /* FIXME unsigned */ int hash ; /* hash value for supernode detection */
1847 int head_column ; /* head of hash bucket */
1848 int first_col ; /* first column in hash bucket */
1849 int tag_mark ; /* marker value for mark array */
1850 int row_mark ; /* Row [row].shared2.mark */
1851 int set_difference ; /* set difference size of row with pivot row */
1852 int min_score ; /* smallest column score */
1853 int col_thickness ; /* "thickness" (no. of columns in a supercol) */
1854 int max_mark ; /* maximum value of tag_mark */
1855 int pivot_col_thickness ; /* number of columns represented by pivot col */
1856 int prev_col ; /* Used by Dlist operations. */
1857 int next_col ; /* Used by Dlist operations. */
1858 int ngarbage ; /* number of garbage collections performed */
1860 int debug_d ; /* debug loop counter */
1861 int debug_step = 0 ; /* debug loop counter */
1863 /* === Initialization and clear mark ================================ */
1865 max_mark = Int_MAX - n_col ;
1866 tag_mark = clear_mark (0, max_mark, n_row, Row) ;
1867 min_score = 0 ;
1868 ngarbage = 0 ;
1869 DEBUG1 ("colamd: Ordering, n_col2=%d\n", n_col2) ;
1871 /* === Order the columns ============================================ */
1873 for (k = 0 ; k < n_col2 ; /* 'k' is incremented below */)
1876 if (!NDEBUG)
1878 if (debug_step % 100 == 0)
1880 DEBUG2 ("\n... Step k: %d out of n_col2: %d\n", k, n_col2) ;
1882 else
1884 DEBUG3 ("\n----------Step k: %d out of n_col2: %d\n", k, n_col2) ;
1886 debug_step++ ;
1887 debug_deg_lists (n_row, n_col, Row, Col, head,
1888 min_score, n_col2-k, max_deg) ;
1889 debug_matrix (n_row, n_col, Row, Col, A) ;
1890 } /* NDEBUG */
1892 /* === Select pivot column, and order it ============================ */
1894 /* make sure degree list isn't empty */
1895 ASSERT (min_score >= 0) ;
1896 ASSERT (min_score <= n_col) ;
1897 ASSERT (head [min_score] >= EMPTY) ;
1899 if (!NDEBUG)
1901 for (debug_d = 0 ; debug_d < min_score ; debug_d++)
1903 ASSERT (head [debug_d] == EMPTY) ;
1905 } /* NDEBUG */
1907 /* get pivot column from head of minimum degree list */
1908 while (head [min_score] == EMPTY && min_score < n_col)
1910 min_score++ ;
1912 pivot_col = head [min_score] ;
1913 ASSERT (pivot_col >= 0 && pivot_col <= n_col) ;
1914 next_col = Col [pivot_col].degree_next() ;
1915 head [min_score] = next_col ;
1916 if (next_col != EMPTY)
1918 Col [next_col].prev(EMPTY) ;
1921 ASSERT (COL_IS_ALIVE (Col, pivot_col)) ;
1923 /* remember score for defrag check */
1924 pivot_col_score = Col [pivot_col].score() ;
1926 /* the pivot column is the kth column in the pivot order */
1927 Col [pivot_col].order(k) ;
1929 /* increment order count by column thickness */
1930 pivot_col_thickness = Col [pivot_col].thickness() ;
1931 k += pivot_col_thickness ;
1932 ASSERT (pivot_col_thickness > 0) ;
1933 DEBUG3 ("Pivot col: %d thick %d\n", pivot_col, pivot_col_thickness) ;
1935 /* === Garbage_collection, if necessary ============================= */
1937 needed_memory = MIN (pivot_col_score, n_col - k) ;
1938 if (pfree + needed_memory >= Alen)
1940 pfree = garbage_collection (n_row, n_col, Row, Col, A, pfree) ;
1941 ngarbage++ ;
1942 /* after garbage collection we will have enough */
1943 ASSERT (pfree + needed_memory < Alen) ;
1944 /* garbage collection has wiped out the Row[].shared2.mark array */
1945 tag_mark = clear_mark (0, max_mark, n_row, Row) ;
1947 if (!NDEBUG)
1949 debug_matrix (n_row, n_col, Row, Col, A) ;
1950 } /* NDEBUG */
1953 /* === Compute pivot row pattern ==================================== */
1955 /* get starting location for this new merged row */
1956 pivot_row_start = pfree ;
1958 /* initialize new row counts to zero */
1959 pivot_row_degree = 0 ;
1961 /* tag pivot column as having been visited so it isn't included */
1962 /* in merged pivot row */
1963 Col [pivot_col].thickness(-pivot_col_thickness) ;
1965 /* pivot row is the union of all rows in the pivot column pattern */
1966 cp = Col [pivot_col].start ;
1967 cp_end = cp + Col [pivot_col].length ;
1968 while (cp < cp_end)
1970 /* get a row */
1971 row = A [cp++] ;
1972 DEBUG4 ("Pivot col pattern %d %d\n", ROW_IS_ALIVE (Row, row) ? 1 : 0, row) ;
1973 /* skip if row is dead */
1974 if (ROW_IS_ALIVE (Row, row))
1976 rp = Row [row].start ;
1977 rp_end = rp + Row [row].length ;
1978 while (rp < rp_end)
1980 /* get a column */
1981 col = A [rp++] ;
1982 /* add the column, if alive and untagged */
1983 col_thickness = Col [col].thickness() ;
1984 if (col_thickness > 0 && COL_IS_ALIVE (Col, col))
1986 /* tag column in pivot row */
1987 Col [col].thickness(-col_thickness) ;
1988 ASSERT (pfree < Alen) ;
1989 /* place column in pivot row */
1990 A [pfree++] = col ;
1991 pivot_row_degree += col_thickness ;
1997 /* clear tag on pivot column */
1998 Col [pivot_col].thickness(pivot_col_thickness) ;
1999 max_deg = MAX (max_deg, pivot_row_degree) ;
2001 if (!NDEBUG)
2003 DEBUG3 ("check2\n") ;
2004 debug_mark (n_row, Row, tag_mark, max_mark) ;
2005 } /* NDEBUG */
2007 /* === Kill all rows used to construct pivot row ==================== */
2009 /* also kill pivot row, temporarily */
2010 cp = Col [pivot_col].start ;
2011 cp_end = cp + Col [pivot_col].length ;
2012 while (cp < cp_end)
2014 /* may be killing an already dead row */
2015 row = A [cp++] ;
2016 DEBUG3 ("Kill row in pivot col: %d\n", row) ;
2017 KILL_ROW (Row, row) ;
2020 /* === Select a row index to use as the new pivot row =============== */
2022 pivot_row_length = pfree - pivot_row_start ;
2023 if (pivot_row_length > 0)
2025 /* pick the "pivot" row arbitrarily (first row in col) */
2026 pivot_row = A [Col [pivot_col].start] ;
2027 DEBUG3 ("Pivotal row is %d\n", pivot_row) ;
2029 else
2031 /* there is no pivot row, since it is of zero length */
2032 pivot_row = EMPTY ;
2033 ASSERT (pivot_row_length == 0) ;
2035 ASSERT (Col [pivot_col].length > 0 || pivot_row_length == 0) ;
2037 /* === Approximate degree computation =============================== */
2039 /* Here begins the computation of the approximate degree. The column */
2040 /* score is the sum of the pivot row "length", plus the size of the */
2041 /* set differences of each row in the column minus the pattern of the */
2042 /* pivot row itself. The column ("thickness") itself is also */
2043 /* excluded from the column score (we thus use an approximate */
2044 /* external degree). */
2046 /* The time taken by the following code (compute set differences, and */
2047 /* add them up) is proportional to the size of the data structure */
2048 /* being scanned - that is, the sum of the sizes of each column in */
2049 /* the pivot row. Thus, the amortized time to compute a column score */
2050 /* is proportional to the size of that column (where size, in this */
2051 /* context, is the column "length", or the number of row indices */
2052 /* in that column). The number of row indices in a column is */
2053 /* monotonically non-decreasing, from the length of the original */
2054 /* column on input to colamd. */
2056 /* === Compute set differences ====================================== */
2058 DEBUG3 ("** Computing set differences phase. **\n") ;
2060 /* pivot row is currently dead - it will be revived later. */
2062 DEBUG3 ("Pivot row: ") ;
2063 /* for each column in pivot row */
2064 rp = pivot_row_start ;
2065 rp_end = rp + pivot_row_length ;
2066 while (rp < rp_end)
2068 col = A [rp++] ;
2069 ASSERT (COL_IS_ALIVE (Col, col) && col != pivot_col) ;
2070 DEBUG3 ("Col: %d\n", col) ;
2072 /* clear tags used to construct pivot row pattern */
2073 col_thickness = -Col [col].thickness() ;
2074 ASSERT (col_thickness > 0) ;
2075 Col [col].thickness(col_thickness) ;
2077 /* === Remove column from degree list =========================== */
2079 cur_score = Col [col].score() ;
2080 prev_col = Col [col].prev() ;
2081 next_col = Col [col].degree_next() ;
2082 ASSERT (cur_score >= 0) ;
2083 ASSERT (cur_score <= n_col) ;
2084 ASSERT (cur_score >= EMPTY) ;
2085 if (prev_col == EMPTY)
2087 head [cur_score] = next_col ;
2089 else
2091 Col [prev_col].degree_next(next_col) ;
2093 if (next_col != EMPTY)
2095 Col [next_col].prev(prev_col) ;
2098 /* === Scan the column ========================================== */
2100 cp = Col [col].start ;
2101 cp_end = cp + Col [col].length ;
2102 while (cp < cp_end)
2104 /* get a row */
2105 row = A [cp++] ;
2106 row_mark = Row [row].mark() ;
2107 /* skip if dead */
2108 if (ROW_IS_MARKED_DEAD (row_mark))
2110 continue ;
2112 ASSERT (row != pivot_row) ;
2113 set_difference = row_mark - tag_mark ;
2114 /* check if the row has been seen yet */
2115 if (set_difference < 0)
2117 ASSERT (Row [row].degree() <= max_deg) ;
2118 set_difference = Row [row].degree() ;
2120 /* subtract column thickness from this row's set difference */
2121 set_difference -= col_thickness ;
2122 ASSERT (set_difference >= 0) ;
2123 /* absorb this row if the set difference becomes zero */
2124 if (set_difference == 0 && aggressive != 0)
2126 DEBUG3 ("aggressive absorption. Row: %d\n", row) ;
2127 KILL_ROW (Row, row) ;
2129 else
2131 /* save the new mark */
2132 Row [row].mark(set_difference + tag_mark) ;
2137 if (!NDEBUG)
2139 debug_deg_lists (n_row, n_col, Row, Col, head,
2140 min_score, n_col2-k-pivot_row_degree, max_deg) ;
2141 } /* NDEBUG */
2143 /* === Add up set differences for each column ======================= */
2145 DEBUG3 ("** Adding set differences phase. **\n") ;
2147 /* for each column in pivot row */
2148 rp = pivot_row_start ;
2149 rp_end = rp + pivot_row_length ;
2150 while (rp < rp_end)
2152 /* get a column */
2153 col = A [rp++] ;
2154 ASSERT (COL_IS_ALIVE (Col, col) && col != pivot_col) ;
2155 hash = 0 ;
2156 cur_score = 0 ;
2157 cp = Col [col].start ;
2158 /* compact the column */
2159 new_cp = cp ;
2160 cp_end = cp + Col [col].length ;
2162 DEBUG4 ("Adding set diffs for Col: %d.\n", col) ;
2164 while (cp < cp_end)
2166 /* get a row */
2167 row = A [cp++] ;
2168 ASSERT(row >= 0 && row < n_row) ;
2169 row_mark = Row [row].mark() ;
2170 /* skip if dead */
2171 if (ROW_IS_MARKED_DEAD (row_mark))
2173 DEBUG4 (" Row %d, dead\n", row) ;
2174 continue ;
2176 DEBUG4 (" Row %d, set diff %d\n", row, row_mark-tag_mark);
2177 ASSERT (row_mark >= tag_mark) ;
2178 /* compact the column */
2179 A [new_cp++] = row ;
2180 /* compute hash function */
2181 hash += row ;
2182 /* add set difference */
2183 cur_score += row_mark - tag_mark ;
2184 /* integer overflow... */
2185 cur_score = MIN (cur_score, n_col) ;
2188 /* recompute the column's length */
2189 Col [col].length = (int) (new_cp - Col [col].start) ;
2191 /* === Further mass elimination ================================= */
2193 if (Col [col].length == 0)
2195 DEBUG4 ("further mass elimination. Col: %d\n", col) ;
2196 /* nothing left but the pivot row in this column */
2197 KILL_PRINCIPAL_COL (Col, col) ;
2198 pivot_row_degree -= Col [col].thickness() ;
2199 ASSERT (pivot_row_degree >= 0) ;
2200 /* order it */
2201 Col [col].order(k) ;
2202 /* increment order count by column thickness */
2203 k += Col [col].thickness() ;
2205 else
2207 /* === Prepare for supercolumn detection ==================== */
2209 DEBUG4 ("Preparing supercol detection for Col: %d.\n", col) ;
2211 /* save score so far */
2212 Col [col].score(cur_score) ;
2214 /* add column to hash table, for supercolumn detection */
2215 hash %= n_col + 1 ;
2217 DEBUG4 (" Hash = %d, n_col = %d.\n", hash, n_col) ;
2218 ASSERT (((int) hash) <= n_col) ;
2220 head_column = head [hash] ;
2221 if (head_column > EMPTY)
2223 /* degree list "hash" is non-empty, use prev (shared3) of */
2224 /* first column in degree list as head of hash bucket */
2225 first_col = Col [head_column].headhash() ;
2226 Col [head_column].headhash(col) ;
2228 else
2230 /* degree list "hash" is empty, use head as hash bucket */
2231 first_col = - (head_column + 2) ;
2232 head [hash] = - (col + 2) ;
2234 Col [col].hash_next(first_col) ;
2236 /* save hash function in Col [col].shared3.hash */
2237 Col [col].hash( (int) hash ) ;
2238 ASSERT (COL_IS_ALIVE (Col, col)) ;
2242 /* The approximate external column degree is now computed. */
2244 /* === Supercolumn detection ======================================== */
2246 DEBUG3 ("** Supercolumn detection phase. **\n") ;
2248 if (!NDEBUG)
2250 detect_super_cols (n_col, Row,
2251 Col, A, head, pivot_row_start, pivot_row_length) ;
2252 } else {
2253 detect_super_cols (Col, A, head, pivot_row_start, pivot_row_length) ;
2256 /* === Kill the pivotal column ====================================== */
2258 KILL_PRINCIPAL_COL (Col, pivot_col) ;
2260 /* === Clear mark =================================================== */
2262 tag_mark = clear_mark (tag_mark+max_deg+1, max_mark, n_row, Row) ;
2264 if (!NDEBUG)
2266 DEBUG3 ("check3\n") ;
2267 debug_mark (n_row, Row, tag_mark, max_mark) ;
2268 } /* NDEBUG */
2270 /* === Finalize the new pivot row, and column scores ================ */
2272 DEBUG3 ("** Finalize scores phase. **\n") ;
2274 /* for each column in pivot row */
2275 rp = pivot_row_start ;
2276 /* compact the pivot row */
2277 new_rp = rp ;
2278 rp_end = rp + pivot_row_length ;
2279 while (rp < rp_end)
2281 col = A [rp++] ;
2282 /* skip dead columns */
2283 if (COL_IS_DEAD (Col, col))
2285 continue ;
2287 A [new_rp++] = col ;
2288 /* add new pivot row to column */
2289 A [Col [col].start + (Col [col].length++)] = pivot_row ;
2291 /* retrieve score so far and add on pivot row's degree. */
2292 /* (we wait until here for this in case the pivot */
2293 /* row's degree was reduced due to mass elimination). */
2294 cur_score = Col [col].score() + pivot_row_degree ;
2296 /* calculate the max possible score as the number of */
2297 /* external columns minus the 'k' value minus the */
2298 /* columns thickness */
2299 max_score = n_col - k - Col [col].thickness() ;
2301 /* make the score the external degree of the union-of-rows */
2302 cur_score -= Col [col].thickness() ;
2304 /* make sure score is less or equal than the max score */
2305 cur_score = MIN (cur_score, max_score) ;
2306 ASSERT (cur_score >= 0) ;
2308 /* store updated score */
2309 Col [col].score(cur_score) ;
2311 /* === Place column back in degree list ========================= */
2313 ASSERT (min_score >= 0) ;
2314 ASSERT (min_score <= n_col) ;
2315 ASSERT (cur_score >= 0) ;
2316 ASSERT (cur_score <= n_col) ;
2317 ASSERT (head [cur_score] >= EMPTY) ;
2318 next_col = head [cur_score] ;
2319 Col [col].degree_next(next_col) ;
2320 Col [col].prev(EMPTY) ;
2321 if (next_col != EMPTY)
2323 Col [next_col].prev(col) ;
2325 head [cur_score] = col ;
2327 /* see if this score is less than current min */
2328 min_score = MIN (min_score, cur_score) ;
2332 if (!NDEBUG)
2334 debug_deg_lists (n_row, n_col, Row, Col, head,
2335 min_score, n_col2-k, max_deg) ;
2336 } /* NDEBUG */
2338 /* === Resurrect the new pivot row ================================== */
2340 if (pivot_row_degree > 0)
2342 /* update pivot row length to reflect any cols that were killed */
2343 /* during super-col detection and mass elimination */
2344 Row [pivot_row].start = pivot_row_start ;
2345 Row [pivot_row].length = (int) (new_rp - pivot_row_start) ;
2346 ASSERT (Row [pivot_row].length > 0) ;
2347 Row [pivot_row].degree(pivot_row_degree) ;
2348 Row [pivot_row].mark(0) ;
2349 /* pivot row is no longer dead */
2351 DEBUG1 ("Resurrect Pivot_row %d deg: %d\n",
2352 pivot_row, pivot_row_degree) ;
2356 /* === All principal columns have now been ordered ================== */
2358 return (ngarbage) ;
2362 /* ========================================================================== */
2363 /* === order_children ======================================================= */
2364 /* ========================================================================== */
2367 * The find_ordering routine has ordered all of the principal columns (the
2368 * representatives of the supercolumns). The non-principal columns have not
2369 * yet been ordered. This routine orders those columns by walking up the
2370 * parent tree (a column is a child of the column which absorbed it). The
2371 * final permutation vector is then placed in p [0 ... n_col-1], with p [0]
2372 * being the first column, and p [n_col-1] being the last. It doesn't look
2373 * like it at first glance, but be assured that this routine takes time linear
2374 * in the number of columns. Although not immediately obvious, the time
2375 * taken by this routine is O (n_col), that is, linear in the number of
2376 * columns. Not user-callable.
2378 * @param n_col number of columns of A
2379 * @param Col of size n_col+1
2380 * @param p p [0 ... n_col-1] is the column permutation
2382 private static void order_children (int n_col, Colamd_Col[] Col, int[] p)
2384 /* === Local variables ============================================== */
2386 int i ; /* loop counter for all columns */
2387 int c ; /* column index */
2388 int parent ; /* index of column's parent */
2389 int order ; /* column's order */
2391 /* === Order each non-principal column ============================== */
2393 for (i = 0 ; i < n_col ; i++)
2395 /* find an un-ordered non-principal column */
2396 ASSERT (COL_IS_DEAD (Col, i)) ;
2397 if (!COL_IS_DEAD_PRINCIPAL (Col, i) && Col [i].order() == EMPTY)
2399 parent = i ;
2400 /* once found, find its principal parent */
2403 parent = Col [parent].parent() ;
2404 } while (!COL_IS_DEAD_PRINCIPAL (Col, parent)) ;
2406 /* now, order all un-ordered non-principal columns along path */
2407 /* to this parent. collapse tree at the same time */
2408 c = i ;
2409 /* get order of parent */
2410 order = Col [parent].order() ;
2414 ASSERT (Col [c].order() == EMPTY) ;
2416 /* order this column */
2417 Col [c].order(order++) ;
2418 /* collaps tree */
2419 Col [c].parent(parent) ;
2421 /* get immediate parent of this column */
2422 c = Col [c].parent() ;
2424 /* continue until we hit an ordered column. There are */
2425 /* guarranteed not to be anymore unordered columns */
2426 /* above an ordered column */
2427 } while (Col [c].order() == EMPTY) ;
2429 /* re-order the super_col parent to largest order for this group */
2430 Col [parent].order(order) ;
2434 /* === Generate the permutation ===================================== */
2436 for (c = 0 ; c < n_col ; c++)
2438 p [Col [c].order()] = c ;
2443 /* ========================================================================== */
2444 /* === detect_super_cols ==================================================== */
2445 /* ========================================================================== */
2448 * Detects supercolumns by finding matches between columns in the hash buckets.
2449 * Check amongst columns in the set A [row_start ... row_start + row_length-1].
2450 * The columns under consideration are currently *not* in the degree lists,
2451 * and have already been placed in the hash buckets.
2453 * The hash bucket for columns whose hash function is equal to h is stored
2454 * as follows:
2456 * if head [h] is >= 0, then head [h] contains a degree list, so:
2458 * head [h] is the first column in degree bucket h.
2459 * Col [head [h]].headhash gives the first column in hash bucket h.
2461 * otherwise, the degree list is empty, and:
2463 * -(head [h] + 2) is the first column in hash bucket h.
2465 * For a column c in a hash bucket, Col [c].shared3.prev is NOT a "previous
2466 * column" pointer. Col [c].shared3.hash is used instead as the hash number
2467 * for that column. The value of Col [c].shared4.hash_next is the next column
2468 * in the same hash bucket.
2470 * Assuming no, or "few" hash collisions, the time taken by this routine is
2471 * linear in the sum of the sizes (lengths) of each column whose score has
2472 * just been computed in the approximate degree computation.
2473 * Not user-callable.
2475 * @param Col of size n_col+1
2476 * @param A row indices of A
2477 * @param head head of degree lists and hash buckets
2478 * @param row_start pointer to set of columns to check
2479 * @param row_length number of columns to check
2481 private static void detect_super_cols (Colamd_Col[] Col, int[] A,
2482 int[] head, int row_start, int row_length)
2484 detect_super_cols (0, null, Col, A, head, row_start, row_length) ;
2488 * debug version
2490 * @param n_col number of columns of A
2491 * @param Row of size n_row+1
2492 * @param Col of size n_col+1
2493 * @param A row indices of A
2494 * @param head head of degree lists and hash buckets
2495 * @param row_start pointer to set of columns to check
2496 * @param row_length number of columns to check
2498 private static void detect_super_cols (int n_col, Colamd_Row[] Row,
2499 Colamd_Col[] Col, int[] A, int[] head, int row_start, int row_length)
2501 /* === Local variables ============================================== */
2503 int hash ; /* hash value for a column */
2504 int rp ; /* pointer to a row */
2505 int c ; /* a column index */
2506 int super_c ; /* column index of the column to absorb into */
2507 int cp1 ; /* column pointer for column super_c */
2508 int cp2 ; /* column pointer for column c */
2509 int length ; /* length of column super_c */
2510 int prev_c ; /* column preceding c in hash bucket */
2511 int i ; /* loop counter */
2512 int rp_end ; /* pointer to the end of the row */
2513 int col ; /* a column index in the row to check */
2514 int head_column ; /* first column in hash bucket or degree list */
2515 int first_col ; /* first column in hash bucket */
2517 /* === Consider each column in the row ============================== */
2519 rp = row_start ;
2520 rp_end = rp + row_length ;
2521 while (rp < rp_end)
2523 col = A [rp++] ;
2524 if (COL_IS_DEAD (Col, col))
2526 continue ;
2529 /* get hash number for this column */
2530 hash = Col [col].hash() ;
2531 ASSERT (hash <= n_col) ;
2533 /* === Get the first column in this hash bucket ===================== */
2535 head_column = head [hash] ;
2536 if (head_column > EMPTY)
2538 first_col = Col [head_column].headhash() ;
2540 else
2542 first_col = - (head_column + 2) ;
2545 /* === Consider each column in the hash bucket ====================== */
2547 for (super_c = first_col ; super_c != EMPTY ;
2548 super_c = Col [super_c].hash_next())
2550 ASSERT (COL_IS_ALIVE (Col, super_c)) ;
2551 ASSERT (Col [super_c].hash() == hash) ;
2552 length = Col [super_c].length ;
2554 /* prev_c is the column preceding column c in the hash bucket */
2555 prev_c = super_c ;
2557 /* === Compare super_c with all columns after it ================ */
2559 for (c = Col [super_c].hash_next() ;
2560 c != EMPTY ; c = Col [c].hash_next())
2562 ASSERT (c != super_c) ;
2563 ASSERT (COL_IS_ALIVE (Col, c)) ;
2564 ASSERT (Col [c].hash() == hash) ;
2566 /* not identical if lengths or scores are different */
2567 if (Col [c].length != length ||
2568 Col [c].score() != Col [super_c].score())
2570 prev_c = c ;
2571 continue ;
2574 /* compare the two columns */
2575 cp1 = Col [super_c].start ;
2576 cp2 = Col [c].start ;
2578 for (i = 0 ; i < length ; i++)
2580 /* the columns are "clean" (no dead rows) */
2581 ASSERT (ROW_IS_ALIVE (Row, A [cp1])) ;
2582 ASSERT (ROW_IS_ALIVE (Row, A [cp2])) ;
2583 /* row indices will same order for both supercols, */
2584 /* no gather scatter nessasary */
2585 if (A [cp1++] != A [cp2++])
2587 break ;
2591 /* the two columns are different if the for-loop "broke" */
2592 if (i != length)
2594 prev_c = c ;
2595 continue ;
2598 /* === Got it! two columns are identical =================== */
2600 ASSERT (Col [c].score() == Col [super_c].score()) ;
2602 Col [super_c].thickness( Col [super_c].thickness() + Col [c].thickness() ) ;
2603 Col [c].parent(super_c) ;
2604 KILL_NON_PRINCIPAL_COL (Col, c) ;
2605 /* order c later, in order_children() */
2606 Col [c].order(EMPTY) ;
2607 /* remove c from hash bucket */
2608 Col [prev_c].hash_next(Col [c].hash_next()) ;
2612 /* === Empty this hash bucket ======================================= */
2614 if (head_column > EMPTY)
2616 /* corresponding degree list "hash" is not empty */
2617 Col [head_column].headhash(EMPTY) ;
2619 else
2621 /* corresponding degree list "hash" is empty */
2622 head [hash] = EMPTY ;
2628 /* ========================================================================== */
2629 /* === garbage_collection =================================================== */
2630 /* ========================================================================== */
2633 * Defragments and compacts columns and rows in the workspace A. Used when
2634 * all avaliable memory has been used while performing row merging. Returns
2635 * the index of the first free position in A, after garbage collection. The
2636 * time taken by this routine is linear in the size of the array A, which is
2637 * itself linear in the number of nonzeros in the input matrix.
2638 * Not user-callable.
2640 * @param n_row number of rows
2641 * @param n_col number of columns
2642 * @param Row row info
2643 * @param Col column info
2644 * @param A A [0 ... Alen-1] holds the matrix
2645 * @param pfree
2646 * @return the new value of pfree
2648 private static int garbage_collection (int n_row, int n_col,
2649 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int pfree)
2651 /* === Local variables ============================================== */
2653 int psrc ; /* source pointer */
2654 int pdest ; /* destination pointer */
2655 int j ; /* counter */
2656 int r ; /* a row index */
2657 int c ; /* a column index */
2658 int length ; /* length of a row or column */
2660 int debug_rows = 0 ;
2661 if (!NDEBUG)
2663 DEBUG2 ("Defrag..\n") ;
2664 for (psrc = 0 ; psrc < pfree ; psrc++) ASSERT (A [psrc] >= 0) ;
2665 debug_rows = 0 ;
2668 /* === Defragment the columns ======================================= */
2670 pdest = 0 ;
2671 for (c = 0 ; c < n_col ; c++)
2673 if (COL_IS_ALIVE (Col, c))
2675 psrc = Col [c].start ;
2677 /* move and compact the column */
2678 ASSERT (pdest <= psrc) ;
2679 Col [c].start = (int) (pdest - 0) ;
2680 length = Col [c].length ;
2681 for (j = 0 ; j < length ; j++)
2683 r = A [psrc++] ;
2684 if (ROW_IS_ALIVE (Row, r))
2686 A [pdest] = r ;
2689 Col [c].length = (int) (pdest - Col [c].start) ;
2693 /* === Prepare to defragment the rows =============================== */
2695 for (r = 0 ; r < n_row ; r++)
2697 if (ROW_IS_DEAD (Row, r) || (Row [r].length == 0))
2699 /* This row is already dead, or is of zero length. Cannot compact
2700 * a row of zero length, so kill it. NOTE: in the current version,
2701 * there are no zero-length live rows. Kill the row (for the first
2702 * time, or again) just to be safe. */
2703 KILL_ROW (Row, r) ;
2705 else
2707 /* save first column index in Row [r].shared2.first_column */
2708 psrc = Row [r].start ;
2709 Row [r].first_column(A [psrc]) ;
2710 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2711 /* flag the start of the row with the one's complement of row */
2712 A [psrc] = ONES_COMPLEMENT (r) ;
2713 if (!NDEBUG)
2715 debug_rows++ ;
2716 } /* NDEBUG */
2720 /* === Defragment the rows ========================================== */
2722 psrc = pdest ;
2723 while (psrc < pfree)
2725 /* find a negative number ... the start of a row */
2726 if (A [psrc++] < 0)
2728 psrc-- ;
2729 /* get the row index */
2730 r = ONES_COMPLEMENT (A [psrc]) ;
2731 ASSERT (r >= 0 && r < n_row) ;
2732 /* restore first column index */
2733 A [psrc] = Row [r].first_column() ;
2734 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2735 ASSERT (Row [r].length > 0) ;
2736 /* move and compact the row */
2737 ASSERT (pdest <= psrc) ;
2738 Row [r].start = (int) (pdest - 0) ;
2739 length = Row [r].length ;
2740 for (j = 0 ; j < length ; j++)
2742 c = A [psrc++] ;
2743 if (COL_IS_ALIVE (Col, c))
2745 A [pdest++] = c ;
2748 Row [r].length = (int) (pdest - Row [r].start) ;
2749 ASSERT (Row [r].length > 0) ;
2750 if (!NDEBUG)
2752 debug_rows-- ;
2753 } /* NDEBUG */
2756 /* ensure we found all the rows */
2757 ASSERT (debug_rows == 0) ;
2759 /* === Return the new value of pfree ================================ */
2761 return ((int) (pdest - 0)) ;
2765 /* ========================================================================== */
2766 /* === clear_mark =========================================================== */
2767 /* ========================================================================== */
2770 * Clears the Row [].shared2.mark array, and returns the new tag_mark.
2772 * @param tag_mark new value of tag_mark
2773 * @param max_mark max allowed value of tag_mark
2774 * @param n_row number of rows in A
2775 * @param Row Row [0 ... n_row-1].shared2.mark is set to zero
2776 * @return the new value for tag_mark
2778 private static int clear_mark (int tag_mark, int max_mark, int n_row,
2779 Colamd_Row[] Row)
2781 int r ;
2783 if (tag_mark <= 0 || tag_mark >= max_mark)
2785 for (r = 0 ; r < n_row ; r++)
2787 if (ROW_IS_ALIVE (Row, r))
2789 Row [r].mark(0) ;
2792 tag_mark = 1 ;
2795 return (tag_mark) ;
2799 /* ========================================================================== */
2800 /* === print_report ========================================================= */
2801 /* ========================================================================== */
2805 * @param method
2806 * @param stats
2808 private static void print_report (String method, int[] stats)
2810 int i1, i2, i3 ;
2812 PRINTF ("\n%s version %d.%d, %s: ", method,
2813 COLAMD_MAIN_VERSION, COLAMD_SUB_VERSION, COLAMD_DATE) ;
2815 if (stats == null)
2817 PRINTF ("No statistics available.\n") ;
2818 return ;
2821 i1 = stats [COLAMD_INFO1] ;
2822 i2 = stats [COLAMD_INFO2] ;
2823 i3 = stats [COLAMD_INFO3] ;
2825 if (stats [COLAMD_STATUS] >= 0)
2827 PRINTF ("OK. ") ;
2829 else
2831 PRINTF ("ERROR. ") ;
2834 switch (stats [COLAMD_STATUS])
2837 case COLAMD_OK_BUT_JUMBLED:
2839 PRINTF("Matrix has unsorted or duplicate row indices.\n") ;
2841 PRINTF("%s: number of duplicate or out-of-order row indices: %d\n",
2842 method, i3) ;
2844 PRINTF("%s: last seen duplicate or out-of-order row index: %d\n",
2845 method, INDEX (i2)) ;
2847 PRINTF("%s: last seen in column: %d",
2848 method, INDEX (i1)) ;
2850 /* no break - fall through to next case instead */
2852 case COLAMD_OK:
2854 PRINTF("\n") ;
2856 PRINTF("%s: number of dense or empty rows ignored: %d\n",
2857 method, stats [COLAMD_DENSE_ROW]) ;
2859 PRINTF("%s: number of dense or empty columns ignored: %d\n",
2860 method, stats [COLAMD_DENSE_COL]) ;
2862 PRINTF("%s: number of garbage collections performed: %d\n",
2863 method, stats [COLAMD_DEFRAG_COUNT]) ;
2864 break ;
2866 case COLAMD_ERROR_A_not_present:
2868 PRINTF("Array A (row indices of matrix) not present.\n") ;
2869 break ;
2871 case COLAMD_ERROR_p_not_present:
2873 PRINTF("Array p (column pointers for matrix) not present.\n") ;
2874 break ;
2876 case COLAMD_ERROR_nrow_negative:
2878 PRINTF("Invalid number of rows (%d).\n", i1) ;
2879 break ;
2881 case COLAMD_ERROR_ncol_negative:
2883 PRINTF("Invalid number of columns (%d).\n", i1) ;
2884 break ;
2886 case COLAMD_ERROR_nnz_negative:
2888 PRINTF("Invalid number of nonzero entries (%d).\n", i1) ;
2889 break ;
2891 case COLAMD_ERROR_p0_nonzero:
2893 PRINTF("Invalid column pointer, p [0] = %d, must be zero.\n", i1);
2894 break ;
2896 case COLAMD_ERROR_A_too_small:
2898 PRINTF("Array A too small.\n") ;
2899 PRINTF(" Need Alen >= %d, but given only Alen = %d.\n",
2900 i1, i2) ;
2901 break ;
2903 case COLAMD_ERROR_col_length_negative:
2905 PRINTF
2906 ("Column %d has a negative number of nonzero entries (%d).\n",
2907 INDEX (i1), i2) ;
2908 break ;
2910 case COLAMD_ERROR_row_index_out_of_bounds:
2912 PRINTF
2913 ("Row index (row %d) out of bounds (%d to %d) in column %d.\n",
2914 INDEX (i2), INDEX (0), INDEX (i3-1), INDEX (i1)) ;
2915 break ;
2917 case COLAMD_ERROR_out_of_memory:
2919 PRINTF("Out of memory.\n") ;
2920 break ;
2922 /* v2.4: internal-error case deleted */
2927 /* ========================================================================== */
2928 /* === colamd debugging routines ============================================ */
2929 /* ========================================================================== */
2932 * At this point, all empty rows and columns are dead. All live columns
2933 * are "clean" (containing no dead rows) and simplicial (no supercolumns
2934 * yet). Rows may contain dead columns, but all live rows contain at
2935 * least one live column.
2937 * @param n_row
2938 * @param n_col
2939 * @param Row
2940 * @param Col
2941 * @param A
2942 * @param n_col2
2944 private static void debug_structures (int n_row, int n_col,
2945 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int n_col2)
2947 if (!NDEBUG)
2949 /* === Local variables ============================================== */
2951 int i ;
2952 int c ;
2953 int cp ;
2954 int cp_end ;
2955 int len ;
2956 int score ;
2957 int r ;
2958 int rp ;
2959 int rp_end ;
2960 int deg ;
2962 /* === Check A, Row, and Col ======================================== */
2964 for (c = 0 ; c < n_col ; c++)
2966 if (COL_IS_ALIVE (Col, c))
2968 len = Col [c].length ;
2969 score = Col [c].score() ;
2970 DEBUG4 ("initial live col %5d %5d %5d\n", c, len, score) ;
2971 ASSERT (len > 0) ;
2972 ASSERT (score >= 0) ;
2973 ASSERT (Col [c].thickness() == 1) ;
2974 cp = Col [c].start ;
2975 cp_end = cp + len ;
2976 while (cp < cp_end)
2978 r = A [cp++] ;
2979 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2982 else
2984 i = Col [c].order() ;
2985 ASSERT (i >= n_col2 && i < n_col) ;
2989 for (r = 0 ; r < n_row ; r++)
2991 if (ROW_IS_ALIVE (Row, r))
2993 i = 0 ;
2994 len = Row [r].length ;
2995 deg = Row [r].degree() ;
2996 ASSERT (len > 0) ;
2997 ASSERT (deg > 0) ;
2998 rp = Row [r].start ;
2999 rp_end = rp + len ;
3000 while (rp < rp_end)
3002 c = A [rp++] ;
3003 if (COL_IS_ALIVE (Col, c))
3005 i++ ;
3008 ASSERT (i > 0) ;
3015 /* ========================================================================== */
3016 /* === debug_deg_lists ====================================================== */
3017 /* ========================================================================== */
3020 * Prints the contents of the degree lists. Counts the number of columns
3021 * in the degree list and compares it to the total it should have. Also
3022 * checks the row degrees.
3024 * @param n_row
3025 * @param n_col
3026 * @param Row
3027 * @param Col
3028 * @param head
3029 * @param min_score
3030 * @param should
3031 * @param max_deg
3033 private static void debug_deg_lists (int n_row, int n_col,
3034 Colamd_Row[] Row, Colamd_Col[] Col, int[] head, int min_score,
3035 int should, int max_deg)
3037 if (!NDEBUG)
3039 /* === Local variables ============================================== */
3041 int deg ;
3042 int col ;
3043 int have ;
3044 int row ;
3046 /* === Check the degree lists ======================================= */
3048 if (n_col > 10000 && colamd_debug <= 0)
3050 return ;
3052 have = 0 ;
3053 DEBUG4 ("Degree lists: %d\n", min_score) ;
3054 for (deg = 0 ; deg <= n_col ; deg++)
3056 col = head [deg] ;
3057 if (col == EMPTY)
3059 continue ;
3061 DEBUG4 ("%d:", deg) ;
3062 while (col != EMPTY)
3064 DEBUG4 (" %d", col) ;
3065 have += Col [col].thickness() ;
3066 ASSERT (COL_IS_ALIVE (Col, col)) ;
3067 col = Col [col].degree_next() ;
3069 DEBUG4 ("\n") ;
3071 DEBUG4 ("should %d have %d\n", should, have) ;
3072 ASSERT (should == have) ;
3074 /* === Check the row degrees ======================================== */
3076 if (n_row > 10000 && colamd_debug <= 0)
3078 return ;
3080 for (row = 0 ; row < n_row ; row++)
3082 if (ROW_IS_ALIVE (Row, row))
3084 ASSERT (Row [row].degree() <= max_deg) ;
3091 /* ========================================================================== */
3092 /* === debug_mark =========================================================== */
3093 /* ========================================================================== */
3096 * Ensures that the tag_mark is less that the maximum and also ensures that
3097 * each entry in the mark array is less than the tag mark.
3099 * @param n_row
3100 * @param Row
3101 * @param tag_mark
3102 * @param max_mark
3104 private static void debug_mark (int n_row, Colamd_Row[] Row, int tag_mark,
3105 int max_mark)
3107 if (!NDEBUG)
3109 /* === Local variables ============================================== */
3111 int r ;
3113 /* === Check the Row marks ========================================== */
3115 ASSERT (tag_mark > 0 && tag_mark <= max_mark) ;
3116 if (n_row > 10000 && colamd_debug <= 0)
3118 return ;
3120 for (r = 0 ; r < n_row ; r++)
3122 ASSERT (Row [r].mark() < tag_mark) ;
3128 /* ========================================================================== */
3129 /* === debug_matrix ========================================================= */
3130 /* ========================================================================== */
3133 * Prints out the contents of the columns and the rows.
3135 * @param n_row
3136 * @param n_col
3137 * @param Row
3138 * @param Col
3139 * @param A
3141 private static void debug_matrix (int n_row, int n_col,
3142 Colamd_Row[] Row, Colamd_Col[] Col, int[] A)
3144 if (!NDEBUG)
3146 /* === Local variables ============================================== */
3148 int r ;
3149 int c ;
3150 int rp ;
3151 int rp_end ;
3152 int cp ;
3153 int cp_end ;
3155 /* === Dump the rows and columns of the matrix ====================== */
3157 if (colamd_debug < 3)
3159 return ;
3161 DEBUG3 ("DUMP MATRIX:\n") ;
3162 for (r = 0 ; r < n_row ; r++)
3164 DEBUG3 ("Row %d alive? %d\n", r, ROW_IS_ALIVE (Row, r) ? 1 : 0) ;
3165 if (ROW_IS_DEAD (Row, r))
3167 continue ;
3169 DEBUG3 ("start %d length %d degree %d\n",
3170 Row [r].start, Row [r].length, Row [r].degree()) ;
3171 rp = Row [r].start ;
3172 rp_end = rp + Row [r].length ;
3173 while (rp < rp_end)
3175 c = A [rp++] ;
3176 DEBUG4 (" %d col %d\n", COL_IS_ALIVE (Col, c) ? 1 : 0, c) ;
3180 for (c = 0 ; c < n_col ; c++)
3182 DEBUG3 ("Col %d alive? %d\n", c, COL_IS_ALIVE (Col, c) ? 1 : 0) ;
3183 if (COL_IS_DEAD (Col, c))
3185 continue ;
3187 DEBUG3 ("start %d length %d shared1 %d shared2 %d\n",
3188 Col [c].start, Col [c].length,
3189 Col [c].thickness(), Col [c].score()) ;
3190 cp = Col [c].start ;
3191 cp_end = cp + Col [c].length ;
3192 while (cp < cp_end)
3194 r = A [cp++] ;
3195 DEBUG4 (" %d row %d\n", ROW_IS_ALIVE (Row, r) ? 1 : 0, r) ;
3201 @SuppressWarnings("unused")
3202 private static void colamd_get_debug (String method)
3204 if (!NDEBUG)
3206 File f ;
3207 colamd_debug = 0 ; /* no debug printing */
3208 f = new File("debug") ;
3209 if (!f.exists())
3211 colamd_debug = 0 ;
3213 else
3215 try {
3216 FileReader fr ;
3217 fr = new FileReader(f) ;
3218 BufferedReader br ;
3219 br = new BufferedReader(fr) ;
3220 colamd_debug = Integer.valueOf( br.readLine() ) ;
3221 br.close() ;
3222 fr.close() ;
3223 } catch (IOException e) {
3224 PRINTF ("%s: AMD_debug_init, " +
3225 "error reading debug.amd file", method) ;
3228 DEBUG0 ("%s: debug version, D = %d (THIS WILL BE SLOW!)\n",
3229 method, colamd_debug) ;