Resolving errors.
[COLAMDJ.git] / src / main / java / edu / ufl / cise / colamd / tdouble / Dcolamd.java
blob2f56078b6cc60091e78a1a5b96ef5788e5616881
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 /* === Row and Column structures ============================================ */
178 /* ========================================================================== */
180 /* User code that makes use of the colamd/symamd routines need not directly */
181 /* reference these structures. They are used only for colamd_recommended. */
183 public class Colamd_Col {
185 /* index for A of first row in this column, or DEAD */
186 /* if column is dead */
187 public int start ;
188 /* number of rows in this column */
189 public int length ;
191 /* number of original columns represented by this */
192 /* col, if the column is alive */
193 public int thickness ;
194 /* parent in parent tree super-column structure, if */
195 /* the column is dead */
196 public int parent ;
198 /* the score used to maintain heap, if col is alive */
199 public int score ;
200 /* pivot ordering of this column, if col is dead */
201 public int order ;
203 /* head of a hash bucket, if col is at the head of */
204 /* a degree list */
205 public int headhash ;
206 /* hash value, if col is not in a degree list */
207 public int hash ;
208 /* previous column in degree list, if col is in a */
209 /* degree list (but not at the head of a degree list) */
210 public int prev ;
212 /* next column, if col is in a degree list */
213 public int degree_next ;
214 /* next column, if col is in a hash list */
215 public int hash_next ;
219 public class Colamd_Row {
221 /* index for A of first col in this row */
222 public int start ;
223 /* number of principal columns in this row */
224 public int length ;
226 /* number of principal & non-principal columns in row */
227 public int degree ;
228 /* used as a row pointer in init_rows_cols () */
229 public int p ;
231 /* for computing set differences and marking dead rows*/
232 public int mark ;
233 /* first column in row (used in garbage collection) */
234 public int first_column ;
238 /* ========================================================================== */
239 /* === Definitions ========================================================== */
240 /* ========================================================================== */
242 protected static double sqrt (double a)
244 return Math.sqrt (a) ;
247 private static final int MAX (int a, int b)
249 return (((a) > (b)) ? (a) : (b)) ;
253 private static final int MIN (int a, int b)
255 return (((a) < (b)) ? (a) : (b)) ;
258 private static final double MAX (double a, double b)
260 return (((a) > (b)) ? (a) : (b)) ;
264 private static int DENSE_DEGREE (double alpha, int n)
266 return ((int) MAX (16.0, (alpha) * sqrt ((double) (n)))) ;
269 private static int ONES_COMPLEMENT (int r)
271 return (-(r)-1) ;
274 private static final int TRUE = (1) ;
275 private static final int FALSE = (0) ;
277 private static final int EMPTY = (-1) ;
279 /* Row and column status */
280 private static final int ALIVE = (0) ;
281 private static final int DEAD = (-1) ;
283 /* Column status */
284 private static final int DEAD_PRINCIPAL = (-1) ;
285 private static final int DEAD_NON_PRINCIPAL = (-2) ;
287 /* Row and column status update and checking. */
288 private static boolean ROW_IS_DEAD(Colamd_Row[] Row, int r)
290 return ROW_IS_MARKED_DEAD (Row [r].mark) ;
292 private static boolean ROW_IS_MARKED_DEAD(int row_mark)
294 return (row_mark < ALIVE) ;
296 private static boolean ROW_IS_ALIVE(Colamd_Row[] Row, int r)
298 return (Row [r].mark >= ALIVE) ;
300 private static boolean COL_IS_DEAD(Colamd_Col[] Col, int c)
302 return (Col [c].start < ALIVE) ;
304 private static boolean COL_IS_ALIVE(Colamd_Col[] Col, int c)
306 return (Col [c].start >= ALIVE) ;
308 private static boolean COL_IS_DEAD_PRINCIPAL(Colamd_Col[] Col, int c)
310 return (Col [c].start == DEAD_PRINCIPAL) ;
312 private static void KILL_ROW(Colamd_Row[] Row, int r)
314 Row [r].mark = DEAD ;
316 private static void KILL_PRINCIPAL_COL(Colamd_Col[] Col, int c)
318 Col [c].start = DEAD_PRINCIPAL ;
320 private static void KILL_NON_PRINCIPAL_COL(Colamd_Col[] Col, int c)
322 Col [c].start = DEAD_NON_PRINCIPAL ;
325 /* ========================================================================== */
326 /* === Colamd reporting mechanism =========================================== */
327 /* ========================================================================== */
330 * In Java, matrices are 0-based and indices are reported as such
331 * in *_report.
333 private static int INDEX (int i)
335 return (i) ;
339 * All output goes through PRINTF.
341 private static void PRINTF (String format, Object... args)
343 if (!NPRINT)
345 System.out.printf (format, args) ;
349 /** debug print level */
350 private static int colamd_debug = 0 ;
352 protected static void DEBUG0 (String format, Object... args)
354 if (!NDEBUG)
356 PRINTF (format, args) ;
360 protected static void DEBUG1(String format, Object... args)
362 if (!NDEBUG)
364 if (colamd_debug >= 1) PRINTF (format, args) ;
368 protected static void DEBUG2(String format, Object... args)
370 if (!NDEBUG)
372 if (colamd_debug >= 2) PRINTF (format, args) ;
376 protected static void DEBUG3(String format, Object... args)
378 if (!NDEBUG)
380 if (colamd_debug >= 3) PRINTF (format, args) ;
384 protected static void DEBUG4(String format, Object... args)
386 if (!NDEBUG)
388 if (colamd_debug >= 4) PRINTF (format, args) ;
392 protected static void ASSERT (boolean a)
394 if (!NDEBUG)
396 assert a ;
400 protected static void ASSERT (int a)
402 ASSERT (a != 0) ;
405 /* ========================================================================== */
406 /* === USER-CALLABLE ROUTINES: ============================================== */
407 /* ========================================================================== */
409 /* ========================================================================== */
410 /* === colamd_recommended =================================================== */
411 /* ========================================================================== */
414 The colamd_recommended routine returns the suggested size for Alen. This
415 value has been determined to provide good balance between the number of
416 garbage collections and the memory requirements for colamd. If any
417 argument is negative, or if integer overflow occurs, a 0 is returned as an
418 error condition. 2*nnz space is required for the row and column
419 indices of the matrix. COLAMD_C (n_col) + COLAMD_R (n_row) space is
420 required for the Col and Row arrays, respectively, which are internal to
421 colamd (roughly 6*n_col + 4*n_row). An additional n_col space is the
422 minimal amount of "elbow room", and nnz/5 more space is recommended for
423 run time efficiency.
425 Alen is approximately 2.2*nnz + 7*n_col + 4*n_row + 10.
427 This function is not needed when using symamd.
431 * add two values of type int, and check for integer overflow
433 private static int t_add (int a, int b, int[] ok)
435 (ok[0]) = (ok[0] != 0) && ((a + b) >= MAX (a,b)) ? 1 : 0;
436 return ((ok[0] != 0) ? (a + b) : 0) ;
440 * compute a*k where k is a small integer, and check for integer overflow
442 static int t_mult (int a, int k, int[] ok)
444 int i, s = 0 ;
445 for (i = 0 ; i < k ; i++)
447 s = t_add (s, a, ok) ;
449 return (s) ;
453 * size of the Col and Row structures
455 private static int COLAMD_C(int n_col, int[] ok)
457 // return ((t_mult (t_add (n_col, 1, ok), sizeof (Colamd_Col), ok) / sizeof (int))) ;
458 return t_add (n_col, 1, ok) ;
461 private static int COLAMD_R(int n_row, int[] ok)
463 // return ((t_mult (t_add (n_row, 1, ok), sizeof (Colamd_Row), ok) / sizeof (int))) ;
464 return t_add (n_row, 1, ok) ;
469 * @param nnz number of nonzeros in A. This must
470 * be the same value as p [n_col] in the call to
471 * colamd - otherwise you will get a wrong value
472 * of the recommended memory to use.
473 * @param n_row number of rows in A
474 * @param n_col number of columns in A
475 * @return recommended value of Alen. Returns 0
476 * if any input argument is negative. The use of this routine
477 * is optional. Not needed for symamd, which dynamically allocates
478 * its own memory.
480 public static int COLAMD_recommended (int nnz, int n_row, int n_col)
482 int s, c, r ;
483 int[] ok = new int [] { TRUE } ;
484 if (nnz < 0 || n_row < 0 || n_col < 0)
486 return (0) ;
488 s = t_mult (nnz, 2, ok) ; /* 2*nnz */
489 // c = COLAMD_C (n_col, ok) ; /* size of column structures */
490 // r = COLAMD_R (n_row, ok) ; /* size of row structures */
491 // s = t_add (s, c, ok) ;
492 // s = t_add (s, r, ok) ;
493 s = t_add (s, n_col, ok) ; /* elbow room */
494 s = t_add (s, nnz/5, ok) ; /* elbow room */
495 ok[0] = (s < Int_MAX) ? 1 : 0;
496 return (ok[0] != 0 ? s : 0) ;
500 /* ========================================================================== */
501 /* === colamd_set_defaults ================================================== */
502 /* ========================================================================== */
505 The colamd_set_defaults routine sets the default values of the user-
506 controllable parameters for colamd and symamd:
508 Colamd: rows with more than max (16, knobs [0] * sqrt (n_col))
509 entries are removed prior to ordering. Columns with more than
510 max (16, knobs [1] * sqrt (MIN (n_row,n_col))) entries are removed
511 prior to ordering, and placed last in the output column ordering.
513 Symamd: Rows and columns with more than max (16, knobs [0] * sqrt (n))
514 entries are removed prior to ordering, and placed last in the
515 output ordering.
517 knobs [0] dense row control
519 knobs [1] dense column control
521 knobs [2] if nonzero, do aggresive absorption
523 knobs [3..19] unused, but future versions might use this
528 * knobs [0] and knobs [1] control dense row and col detection:
530 * Colamd: rows with more than
531 * max (16, knobs [COLAMD_DENSE_ROW] * sqrt (n_col))
532 * entries are removed prior to ordering. Columns with more than
533 * max (16, knobs [COLAMD_DENSE_COL] * sqrt (MIN (n_row,n_col)))
534 * entries are removed prior to
535 * ordering, and placed last in the output column ordering.
537 * Symamd: uses only knobs [COLAMD_DENSE_ROW], which is knobs [0].
538 * Rows and columns with more than
539 * max (16, knobs [COLAMD_DENSE_ROW] * sqrt (n))
540 * entries are removed prior to ordering, and placed last in the
541 * output ordering.
543 * COLAMD_DENSE_ROW and COLAMD_DENSE_COL are defined as 0 and 1,
544 * respectively, in colamd.h. Default values of these two knobs
545 * are both 10. Currently, only knobs [0] and knobs [1] are
546 * used, but future versions may use more knobs. If so, they will
547 * be properly set to their defaults by the future version of
548 * colamd_set_defaults, so that the code that calls colamd will
549 * not need to change, assuming that you either use
550 * colamd_set_defaults, or pass a (double *) NULL pointer as the
551 * knobs array to colamd or symamd.
553 * knobs [2]: aggressive absorption
555 * knobs [COLAMD_AGGRESSIVE] controls whether or not to do
556 * aggressive absorption during the ordering. Default is TRUE.
558 * @param knobs knob array
560 public static void COLAMD_set_defaults (double[] knobs)
562 /* === Local variables ============================================== */
564 int i ;
566 if (knobs == null || knobs.length == 0)
568 return ; /* no knobs to initialize */
570 for (i = 0 ; i < COLAMD_KNOBS ; i++)
572 knobs [i] = 0 ;
574 knobs [COLAMD_DENSE_ROW] = 10 ;
575 knobs [COLAMD_DENSE_COL] = 10 ;
576 knobs [COLAMD_AGGRESSIVE] = TRUE ; /* default: do aggressive absorption*/
579 /* ========================================================================== */
580 /* === symamd =============================================================== */
581 /* ========================================================================== */
584 * The symamd routine computes an ordering P of a symmetric sparse
585 * matrix A such that the Cholesky factorization PAP' = LL' remains
586 * sparse. It is based on a column ordering of a matrix M constructed
587 * so that the nonzero pattern of M'M is the same as A. The matrix A
588 * is assumed to be symmetric; only the strictly lower triangular part
589 * is accessed. You must pass your selected memory allocator (usually
590 * calloc/free or mxCalloc/mxFree) to symamd, for it to allocate
591 * memory for the temporary matrix M.
593 * @param n number of rows and columns of A. Restriction: n >= 0.
594 * Symamd returns FALSE if n is negative.
595 * @param A an integer array of size nnz, where nnz = p [n].
597 * The row indices of the entries in column c of the matrix are
598 * held in A [(p [c]) ... (p [c+1]-1)]. The row indices in a
599 * given column c need not be in ascending order, and duplicate
600 * row indices may be present. However, symamd will run faster
601 * if the columns are in sorted order with no duplicate entries.
603 * The matrix is 0-based. That is, rows are in the range 0 to
604 * n-1, and columns are in the range 0 to n-1. Symamd
605 * returns FALSE if any row index is out of range.
607 * The contents of A are not modified.
608 * @param an integer array of size n+1. On input, it holds the
609 * "pointers" for the column form of the matrix A. Column c of
610 * the matrix A is held in A [(p [c]) ... (p [c+1]-1)]. The first
611 * entry, p [0], must be zero, and p [c] <= p [c+1] must hold
612 * for all c in the range 0 to n-1. The value p [n] is
613 * thus the total number of entries in the pattern of the matrix A.
614 * Symamd returns FALSE if these conditions are not met.
616 * The contents of p are not modified.
617 * @param perm On output, if symamd returns TRUE, the array perm holds the
618 * permutation P, where perm [0] is the first index in the new
619 * ordering, and perm [n-1] is the last. That is, perm [k] = j
620 * means that row and column j of A is the kth column in PAP',
621 * where k is in the range 0 to n-1 (perm [0] = j means
622 * that row and column j of A are the first row and column in
623 * PAP'). The array is used as a workspace during the ordering,
624 * which is why it must be of length n+1, not just n.
625 * @param knobs parameters (uses defaults if NULL)
626 * @param stats Statistics on the ordering, and error status.
627 * Symamd returns FALSE if stats is not present.
629 * stats [0]: number of dense or empty row and columns ignored
630 * (and ordered last in the output permutation
631 * perm). Note that a row/column can become
632 * "empty" if it contains only "dense" and/or
633 * "empty" columns/rows.
635 * stats [1]: (same as stats [0])
637 * stats [2]: number of garbage collections performed.
639 * stats [3]: status code. < 0 is an error code.
640 * > 1 is a warning or notice.
642 * 0 OK. Each column of the input matrix contained
643 * row indices in increasing order, with no
644 * duplicates.
646 * 1 OK, but columns of input matrix were jumbled
647 * (unsorted columns or duplicate entries). Symamd
648 * had to do some extra work to sort the matrix
649 * first and remove duplicate entries, but it
650 * still was able to return a valid permutation
651 * (return value of symamd was TRUE).
653 * stats [4]: highest numbered column that
654 * is unsorted or has duplicate
655 * entries.
656 * stats [5]: last seen duplicate or
657 * unsorted row index.
658 * stats [6]: number of duplicate or
659 * unsorted row indices.
661 * -1 A is a null pointer
663 * -2 p is a null pointer
665 * -3 (unused, see colamd.c)
667 * -4 n is negative
669 * stats [4]: n
671 * -5 number of nonzeros in matrix is negative
673 * stats [4]: # of nonzeros (p [n]).
675 * -6 p [0] is nonzero
677 * stats [4]: p [0]
679 * -7 (unused)
681 * -8 a column has a negative number of entries
683 * stats [4]: column with < 0 entries
684 * stats [5]: number of entries in col
686 * -9 a row index is out of bounds
688 * stats [4]: column with bad row index
689 * stats [5]: bad row index
690 * stats [6]: n_row, # of rows of matrx
692 * -10 out of memory (unable to allocate temporary
693 * workspace for M or count arrays using the
694 * "allocate" routine passed into symamd).
696 * Future versions may return more statistics in the stats array.
697 * @param allocate pointer to calloc
698 * @param release pointer to free
699 * @return TRUE if OK, FALSE otherwise
701 public static int symamd (int n, int[] A, int[] p, int[] perm,
702 double[] knobs, int[] stats)
704 /* === Local variables ============================================== */
706 int[] count ; /* length of each column of M, and col pointer*/
707 int[] mark ; /* mark array for finding duplicate entries */
708 int[] M ; /* row indices of matrix M */
709 int Mlen ; /* length of M */
710 int n_row ; /* number of rows in M */
711 int nnz ; /* number of entries in A */
712 int i ; /* row index of A */
713 int j ; /* column index of A */
714 int k ; /* row index of M */
715 int mnz ; /* number of nonzeros in M */
716 int pp ; /* index into a column of A */
717 int last_row ; /* last row seen in the current column */
718 int length ; /* number of nonzeros in a column */
720 double[] cknobs = new double[COLAMD_KNOBS] ; /* knobs for colamd */
721 double[] default_knobs = new double[COLAMD_KNOBS] ; /* default knobs for colamd */
723 if (!NDEBUG)
725 // colamd_get_debug ("symamd") ;
728 /* === Check the input arguments ==================================== */
730 if (stats == null)
732 DEBUG0 ("symamd: stats not present\n") ;
733 return (FALSE) ;
735 for (i = 0 ; i < COLAMD_STATS ; i++)
737 stats [i] = 0 ;
739 stats [COLAMD_STATUS] = COLAMD_OK ;
740 stats [COLAMD_INFO1] = -1 ;
741 stats [COLAMD_INFO2] = -1 ;
743 if (A == null)
745 stats [COLAMD_STATUS] = COLAMD_ERROR_A_not_present ;
746 DEBUG0 ("symamd: A not present\n") ;
747 return (FALSE) ;
750 if (p == null) /* p is not present */
752 stats [COLAMD_STATUS] = COLAMD_ERROR_p_not_present ;
753 DEBUG0 ("symamd: p not present\n") ;
754 return (FALSE) ;
757 if (n < 0) /* n must be >= 0 */
759 stats [COLAMD_STATUS] = COLAMD_ERROR_ncol_negative ;
760 stats [COLAMD_INFO1] = n ;
761 DEBUG0 ("symamd: n negative %d\n", n) ;
762 return (FALSE) ;
765 nnz = p [n] ;
766 if (nnz < 0) /* nnz must be >= 0 */
768 stats [COLAMD_STATUS] = COLAMD_ERROR_nnz_negative ;
769 stats [COLAMD_INFO1] = nnz ;
770 DEBUG0 ("symamd: number of entries negative %d\n", nnz) ;
771 return (FALSE) ;
774 if (p [0] != 0)
776 stats [COLAMD_STATUS] = COLAMD_ERROR_p0_nonzero ;
777 stats [COLAMD_INFO1] = p [0] ;
778 DEBUG0 ("symamd: p[0] not zero %d\n", p [0]) ;
779 return (FALSE) ;
782 /* === If no knobs, set default knobs =============================== */
784 if (knobs == null)
786 COLAMD_set_defaults (default_knobs) ;
787 knobs = default_knobs ;
790 /* === Allocate count and mark ====================================== */
794 count = new int [n+1] ;
795 } catch (OutOfMemoryError e) {
796 stats [COLAMD_STATUS] = COLAMD_ERROR_out_of_memory ;
797 DEBUG0 ("symamd: allocate count (size %d) failed\n", n+1) ;
798 return (FALSE) ;
803 mark = new int[n+1] ;
804 } catch (OutOfMemoryError e) {
805 stats [COLAMD_STATUS] = COLAMD_ERROR_out_of_memory ;
806 count = null ;
807 DEBUG0 ("symamd: allocate mark (size %d) failed\n", n+1) ;
808 return (FALSE) ;
811 /* === Compute column counts of M, check if A is valid ============== */
813 stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
815 for (i = 0 ; i < n ; i++)
817 mark [i] = -1 ;
820 for (j = 0 ; j < n ; j++)
822 last_row = -1 ;
824 length = p [j+1] - p [j] ;
825 if (length < 0)
827 /* column pointers must be non-decreasing */
828 stats [COLAMD_STATUS] = COLAMD_ERROR_col_length_negative ;
829 stats [COLAMD_INFO1] = j ;
830 stats [COLAMD_INFO2] = length ;
831 count = null ;
832 mark = null ;
833 DEBUG0 ("symamd: col %d negative length %d\n", j, length) ;
834 return (FALSE) ;
837 for (pp = p [j] ; pp < p [j+1] ; pp++)
839 i = A [pp] ;
840 if (i < 0 || i >= n)
842 /* row index i, in column j, is out of bounds */
843 stats [COLAMD_STATUS] = COLAMD_ERROR_row_index_out_of_bounds ;
844 stats [COLAMD_INFO1] = j ;
845 stats [COLAMD_INFO2] = i ;
846 stats [COLAMD_INFO3] = n ;
847 count = null ;
848 mark = null ;
849 DEBUG0 ("symamd: row %d col %d out of bounds\n", i, j) ;
850 return (FALSE) ;
853 if (i <= last_row || mark [i] == j)
855 /* row index is unsorted or repeated (or both), thus col */
856 /* is jumbled. This is a notice, not an error condition. */
857 stats [COLAMD_STATUS] = COLAMD_OK_BUT_JUMBLED ;
858 stats [COLAMD_INFO1] = j ;
859 stats [COLAMD_INFO2] = i ;
860 (stats [COLAMD_INFO3]) ++ ;
861 DEBUG1 ("symamd: row %d col %d unsorted/duplicate\n", i, j) ;
864 if (i > j && mark [i] != j)
866 /* row k of M will contain column indices i and j */
867 count [i]++ ;
868 count [j]++ ;
871 /* mark the row as having been seen in this column */
872 mark [i] = j ;
874 last_row = i ;
878 /* === Compute column pointers of M ================================= */
880 /* use output permutation, perm, for column pointers of M */
881 perm [0] = 0 ;
882 for (j = 1 ; j <= n ; j++)
884 perm [j] = perm [j-1] + count [j-1] ;
886 for (j = 0 ; j < n ; j++)
888 count [j] = perm [j] ;
891 /* === Construct M ================================================== */
893 mnz = perm [n] ;
894 n_row = mnz / 2 ;
895 Mlen = COLAMD_recommended (mnz, n_row, n) ;
898 M = new int [Mlen] ;
899 DEBUG0 ("symamd: M is %d-by-%d with %d entries, Mlen = %g\n",
900 n_row, n, mnz, Mlen) ;
901 } catch (OutOfMemoryError e) {
902 stats [COLAMD_STATUS] = COLAMD_ERROR_out_of_memory ;
903 count = null ;
904 mark = null;
905 DEBUG0 ("symamd: allocate M (size %g) failed\n", (double) Mlen) ;
906 return (FALSE) ;
909 k = 0 ;
911 if (stats [COLAMD_STATUS] == COLAMD_OK)
913 /* Matrix is OK */
914 for (j = 0 ; j < n ; j++)
916 ASSERT (p [j+1] - p [j] >= 0) ;
917 for (pp = p [j] ; pp < p [j+1] ; pp++)
919 i = A [pp] ;
920 ASSERT (i >= 0 && i < n) ;
921 if (i > j)
923 /* row k of M contains column indices i and j */
924 M [count [i]++] = k ;
925 M [count [j]++] = k ;
926 k++ ;
931 else
933 /* Matrix is jumbled. Do not add duplicates to M. Unsorted cols OK. */
934 DEBUG0 ("symamd: Duplicates in A.\n") ;
935 for (i = 0 ; i < n ; i++)
937 mark [i] = -1 ;
939 for (j = 0 ; j < n ; j++)
941 ASSERT (p [j+1] - p [j] >= 0) ;
942 for (pp = p [j] ; pp < p [j+1] ; pp++)
944 i = A [pp] ;
945 ASSERT (i >= 0 && i < n) ;
946 if (i > j && mark [i] != j)
948 /* row k of M contains column indices i and j */
949 M [count [i]++] = k ;
950 M [count [j]++] = k ;
951 k++ ;
952 mark [i] = j ;
958 /* count and mark no longer needed */
959 count = null ;
960 mark = null ;
961 ASSERT (k == n_row) ;
963 /* === Adjust the knobs for M ======================================= */
965 for (i = 0 ; i < COLAMD_KNOBS ; i++)
967 cknobs [i] = knobs [i] ;
970 /* there are no dense rows in M */
971 cknobs [COLAMD_DENSE_ROW] = -1 ;
972 cknobs [COLAMD_DENSE_COL] = knobs [COLAMD_DENSE_ROW] ;
974 /* === Order the columns of M ======================================= */
976 colamd (n_row, n, Mlen, M, perm, cknobs, stats) ;
978 /* Note that the output permutation is now in perm */
980 /* === get the statistics for symamd from colamd ==================== */
982 /* a dense column in colamd means a dense row and col in symamd */
983 stats [COLAMD_DENSE_ROW] = stats [COLAMD_DENSE_COL] ;
985 /* === Free M ======================================================= */
987 M = null ;
988 DEBUG0 ("symamd: done.\n") ;
989 return (TRUE) ;
992 /* ========================================================================== */
993 /* === colamd =============================================================== */
994 /* ========================================================================== */
997 The colamd routine computes a column ordering Q of a sparse matrix
998 A such that the LU factorization P(AQ) = LU remains sparse, where P is
999 selected via partial pivoting. The routine can also be viewed as
1000 providing a permutation Q such that the Cholesky factorization
1001 (AQ)'(AQ) = LL' remains sparse.
1005 * Computes a column ordering (Q) of A such that P(AQ)=LU or
1006 * (AQ)'AQ=LL' have less fill-in and require fewer floating point
1007 * operations than factorizing the unpermuted matrix A or A'A,
1008 * respectively.
1010 * @param n_row number of rows in A. Restriction: n_row >= 0.
1011 * Colamd returns FALSE if n_row is negative.
1012 * @param n_col number of columns in A. Restriction: n_col >= 0.
1013 * Colamd returns FALSE if n_col is negative.
1014 * @param Alen length of A. Restriction (see note):
1015 * Alen >= 2*nnz + 6*(n_col+1) + 4*(n_row+1) + n_col
1016 * Colamd returns FALSE if these conditions are not met.
1018 * Note: this restriction makes an modest assumption regarding
1019 * the size of the two typedef's structures in colamd.h.
1020 * We do, however, guarantee that
1022 * Alen >= colamd_recommended (nnz, n_row, n_col)
1024 * will be sufficient. Note: the macro version does not check
1025 * for integer overflow, and thus is not recommended. Use
1026 * the colamd_recommended routine instead.
1027 * @param A row indices of A.
1029 * A is an integer array of size Alen. Alen must be at least as
1030 * large as the bare minimum value given above, but this is very
1031 * low, and can result in excessive run time. For best
1032 * performance, we recommend that Alen be greater than or equal to
1033 * colamd_recommended (nnz, n_row, n_col), which adds
1034 * nnz/5 to the bare minimum value given above.
1036 * On input, the row indices of the entries in column c of the
1037 * matrix are held in A [(p [c]) ... (p [c+1]-1)]. The row indices
1038 * in a given column c need not be in ascending order, and
1039 * duplicate row indices may be be present. However, colamd will
1040 * work a little faster if both of these conditions are met
1041 * (Colamd puts the matrix into this format, if it finds that the
1042 * the conditions are not met).
1044 * The matrix is 0-based. That is, rows are in the range 0 to
1045 * n_row-1, and columns are in the range 0 to n_col-1. Colamd
1046 * returns FALSE if any row index is out of range.
1048 * The contents of A are modified during ordering, and are
1049 * undefined on output.
1050 * @param p pointers to columns in A.
1052 * p is an integer array of size n_col+1. On input, it holds the
1053 * "pointers" for the column form of the matrix A. Column c of
1054 * the matrix A is held in A [(p [c]) ... (p [c+1]-1)]. The first
1055 * entry, p [0], must be zero, and p [c] <= p [c+1] must hold
1056 * for all c in the range 0 to n_col-1. The value p [n_col] is
1057 * thus the total number of entries in the pattern of the matrix A.
1058 * Colamd returns FALSE if these conditions are not met.
1060 * On output, if colamd returns TRUE, the array p holds the column
1061 * permutation (Q, for P(AQ)=LU or (AQ)'(AQ)=LL'), where p [0] is
1062 * the first column index in the new ordering, and p [n_col-1] is
1063 * the last. That is, p [k] = j means that column j of A is the
1064 * kth pivot column, in AQ, where k is in the range 0 to n_col-1
1065 * (p [0] = j means that column j of A is the first column in AQ).
1067 * If colamd returns FALSE, then no permutation is returned, and
1068 * p is undefined on output.
1069 * @param knobs parameters (uses defaults if NULL)
1070 * @param stats output statistics and error codes.
1072 * Statistics on the ordering, and error status.
1073 * Colamd returns FALSE if stats is not present.
1075 * stats [0]: number of dense or empty rows ignored.
1077 * stats [1]: number of dense or empty columns ignored (and
1078 * ordered last in the output permutation p)
1079 * Note that a row can become "empty" if it
1080 * contains only "dense" and/or "empty" columns,
1081 * and similarly a column can become "empty" if it
1082 * only contains "dense" and/or "empty" rows.
1084 * stats [2]: number of garbage collections performed.
1085 * This can be excessively high if Alen is close
1086 * to the minimum required value.
1088 * stats [3]: status code. < 0 is an error code.
1089 * > 1 is a warning or notice.
1091 * 0 OK. Each column of the input matrix contained
1092 * row indices in increasing order, with no
1093 * duplicates.
1095 * 1 OK, but columns of input matrix were jumbled
1096 * (unsorted columns or duplicate entries). Colamd
1097 * had to do some extra work to sort the matrix
1098 * first and remove duplicate entries, but it
1099 * still was able to return a valid permutation
1100 * (return value of colamd was TRUE).
1102 * stats [4]: highest numbered column that
1103 * is unsorted or has duplicate
1104 * entries.
1105 * stats [5]: last seen duplicate or
1106 * unsorted row index.
1107 * stats [6]: number of duplicate or
1108 * unsorted row indices.
1110 * -1 A is a null pointer
1112 * -2 p is a null pointer
1114 * -3 n_row is negative
1116 * stats [4]: n_row
1118 * -4 n_col is negative
1120 * stats [4]: n_col
1122 * -5 number of nonzeros in matrix is negative
1124 * stats [4]: number of nonzeros, p [n_col]
1126 * -6 p [0] is nonzero
1128 * stats [4]: p [0]
1130 * -7 A is too small
1132 * stats [4]: required size
1133 * stats [5]: actual size (Alen)
1135 * -8 a column has a negative number of entries
1137 * stats [4]: column with < 0 entries
1138 * stats [5]: number of entries in col
1140 * -9 a row index is out of bounds
1142 * stats [4]: column with bad row index
1143 * stats [5]: bad row index
1144 * stats [6]: n_row, # of rows of matrx
1146 * -10 (unused; see symamd.c)
1148 * -999 (unused; see symamd.c)
1150 * Future versions may return more statistics in the stats array.
1151 * @return TRUE if successful, FALSE otherwise
1153 public static int colamd (int n_row, int n_col, int Alen, int[] A,
1154 int[] p, double[] knobs, int[] stats)
1156 /* === Local variables ============================================== */
1158 int i ; /* loop index */
1159 int nnz ; /* nonzeros in A */
1160 int Row_size ; /* size of Row [], in integers */
1161 int Col_size ; /* size of Col [], in integers */
1162 int need ; /* minimum required length of A */
1163 Colamd_Row[] Row ; /* pointer into A of Row [0..n_row] array */
1164 Colamd_Col[] Col ; /* pointer into A of Col [0..n_col] array */
1165 int[] n_col2 = new int [1] ; /* number of non-dense, non-empty columns */
1166 int[] n_row2 = new int [1] ; /* number of non-dense, non-empty rows */
1167 int ngarbage ; /* number of garbage collections performed */
1168 int[] max_deg = new int [1] ; /* maximum row degree */
1169 double[] default_knobs = new double[COLAMD_KNOBS] ; /* default knobs array */
1170 int aggressive ; /* do aggressive absorption */
1171 int[] ok ;
1173 if (!NDEBUG)
1175 // colamd_get_debug ("colamd") ;
1178 /* === Check the input arguments ==================================== */
1180 if (stats == null)
1182 DEBUG0 ("colamd: stats not present\n") ;
1183 return (FALSE) ;
1185 for (i = 0 ; i < COLAMD_STATS ; i++)
1187 stats [i] = 0 ;
1189 stats [COLAMD_STATUS] = COLAMD_OK ;
1190 stats [COLAMD_INFO1] = -1 ;
1191 stats [COLAMD_INFO2] = -1 ;
1193 if (A == null) /* A is not present */
1195 stats [COLAMD_STATUS] = COLAMD_ERROR_A_not_present ;
1196 DEBUG0 ("colamd: A not present\n") ;
1197 return (FALSE) ;
1200 if (p == null) /* p is not present */
1202 stats [COLAMD_STATUS] = COLAMD_ERROR_p_not_present ;
1203 DEBUG0 ("colamd: p not present\n") ;
1204 return (FALSE) ;
1207 if (n_row < 0) /* n_row must be >= 0 */
1209 stats [COLAMD_STATUS] = COLAMD_ERROR_nrow_negative ;
1210 stats [COLAMD_INFO1] = n_row ;
1211 DEBUG0 ("colamd: nrow negative %d\n", n_row) ;
1212 return (FALSE) ;
1215 if (n_col < 0) /* n_col must be >= 0 */
1217 stats [COLAMD_STATUS] = COLAMD_ERROR_ncol_negative ;
1218 stats [COLAMD_INFO1] = n_col ;
1219 DEBUG0 ("colamd: ncol negative %d\n", n_col) ;
1220 return (FALSE) ;
1223 nnz = p [n_col] ;
1224 if (nnz < 0) /* nnz must be >= 0 */
1226 stats [COLAMD_STATUS] = COLAMD_ERROR_nnz_negative ;
1227 stats [COLAMD_INFO1] = nnz ;
1228 DEBUG0 ("colamd: number of entries negative %d\n", nnz) ;
1229 return (FALSE) ;
1232 if (p [0] != 0)
1234 stats [COLAMD_STATUS] = COLAMD_ERROR_p0_nonzero ;
1235 stats [COLAMD_INFO1] = p [0] ;
1236 DEBUG0 ("colamd: p[0] not zero %d\n", p [0]) ;
1237 return (FALSE) ;
1240 /* === If no knobs, set default knobs =============================== */
1242 if (knobs == null)
1244 COLAMD_set_defaults (default_knobs) ;
1245 knobs = default_knobs ;
1248 aggressive = (knobs [COLAMD_AGGRESSIVE] != FALSE) ? 1 : 0;
1250 /* === Allocate the Row and Col arrays from array A ================= */
1252 ok = new int [] { TRUE } ;
1253 Col_size = COLAMD_C (n_col, ok) ; /* size of Col array of structs */
1254 Row_size = COLAMD_R (n_row, ok) ; /* size of Row array of structs */
1256 /* need = 2*nnz + n_col + Col_size + Row_size ; */
1257 need = t_mult (nnz, 2, ok) ;
1258 need = t_add (need, n_col, ok) ;
1259 // need = t_add (need, Col_size, ok) ;
1260 // need = t_add (need, Row_size, ok) ;
1262 if ((ok[0] == 0) || need > (int) Alen || need > Int_MAX)
1264 /* not enough space in array A to perform the ordering */
1265 stats [COLAMD_STATUS] = COLAMD_ERROR_A_too_small ;
1266 stats [COLAMD_INFO1] = need ;
1267 stats [COLAMD_INFO2] = Alen ;
1268 DEBUG0 ("colamd: Need Alen >= %d, given only Alen = %d\n", need,Alen);
1269 return (FALSE) ;
1272 // Alen -= Col_size + Row_size ;
1273 Col = new Colamd_Col [Col_size] ; //A [Alen] ;
1274 Row = new Colamd_Row [Row_size] ; //A [Alen + Col_size] ;
1276 /* === Construct the row and column data structures ================= */
1278 if (init_rows_cols (n_row, n_col, Row, Col, A, p, stats) == 0)
1280 /* input matrix is invalid */
1281 DEBUG0 ("colamd: Matrix invalid\n") ;
1282 return (FALSE) ;
1285 /* === Initialize scores, kill dense rows/columns =================== */
1287 init_scoring (n_row, n_col, Row, Col, A, p, knobs,
1288 n_row2, n_col2, max_deg) ;
1290 /* === Order the supercolumns ======================================= */
1292 ngarbage = find_ordering (n_row, n_col, Alen, Row, Col, A, p,
1293 n_col2[0], max_deg[0], 2*nnz, aggressive) ;
1295 /* === Order the non-principal columns ============================== */
1297 order_children (n_col, Col, p) ;
1299 /* === Return statistics in stats =================================== */
1301 stats [COLAMD_DENSE_ROW] = n_row - n_row2[0] ;
1302 stats [COLAMD_DENSE_COL] = n_col - n_col2[0] ;
1303 stats [COLAMD_DEFRAG_COUNT] = ngarbage ;
1304 DEBUG0 ("colamd: done.\n") ;
1305 return (TRUE) ;
1309 /* ========================================================================== */
1310 /* === colamd_report ======================================================== */
1311 /* ========================================================================== */
1313 public static void COLAMD_report (int[] stats)
1315 print_report ("colamd", stats) ;
1319 /* ========================================================================== */
1320 /* === symamd_report ======================================================== */
1321 /* ========================================================================== */
1323 public static void SYMAMD_report (int[] stats)
1325 print_report ("symamd", stats) ;
1330 /* ========================================================================== */
1331 /* === NON-USER-CALLABLE ROUTINES: ========================================== */
1332 /* ========================================================================== */
1334 /* There are no user-callable routines beyond this point in the file */
1337 /* ========================================================================== */
1338 /* === init_rows_cols ======================================================= */
1339 /* ========================================================================== */
1342 * Takes the column form of the matrix in A and creates the row form of the
1343 * matrix. Also, row and column attributes are stored in the Col and Row
1344 * structs. If the columns are un-sorted or contain duplicate row indices,
1345 * this routine will also sort and remove duplicate row indices from the
1346 * column form of the matrix. Returns FALSE if the matrix is invalid,
1347 * TRUE otherwise. Not user-callable.
1349 * @param n_row number of rows of A
1350 * @param n_col number of columns of A
1351 * @param Row of size n_row+1
1352 * @param Col of size n_col+1
1353 * @param A row indices of A, of size Alen
1354 * @param p pointers to columns in A, of size n_col+1
1355 * @param stats colamd statistics
1356 * @return TRUE if OK, or FALSE otherwise
1358 private static int init_rows_cols (int n_row, int n_col, Colamd_Row[] Row,
1359 Colamd_Col[] Col, int[] A, int[] p, int[] stats)
1361 /* === Local variables ============================================== */
1363 int col ; /* a column index */
1364 int row ; /* a row index */
1365 int cp ; /* a column pointer */
1366 int cp_end ; /* a pointer to the end of a column */
1367 int rp ; /* a row pointer */
1368 int rp_end ; /* a pointer to the end of a row */
1369 int last_row ; /* previous row */
1371 /* === Initialize columns, and check column pointers ================ */
1373 for (col = 0 ; col < n_col ; col++)
1375 Col [col].start = p [col] ;
1376 Col [col].length = p [col+1] - p [col] ;
1378 if (Col [col].length < 0)
1380 /* column pointers must be non-decreasing */
1381 stats [COLAMD_STATUS] = COLAMD_ERROR_col_length_negative ;
1382 stats [COLAMD_INFO1] = col ;
1383 stats [COLAMD_INFO2] = Col [col].length ;
1384 DEBUG0 ("colamd: col %d length %d < 0\n", col, Col [col].length) ;
1385 return (FALSE) ;
1388 Col [col].thickness = 1 ;
1389 Col [col].score = 0 ;
1390 Col [col].prev = EMPTY ;
1391 Col [col].degree_next = EMPTY ;
1394 /* p [0..n_col] no longer needed, used as "head" in subsequent routines */
1396 /* === Scan columns, compute row degrees, and check row indices ===== */
1398 stats [COLAMD_INFO3] = 0 ; /* number of duplicate or unsorted row indices*/
1400 for (row = 0 ; row < n_row ; row++)
1402 Row [row].length = 0 ;
1403 Row [row].mark = -1 ;
1406 for (col = 0 ; col < n_col ; col++)
1408 last_row = -1 ;
1410 cp = A [p [col]] ;
1411 cp_end = A [p [col+1]] ;
1413 while (cp < cp_end)
1415 row = cp++ ;
1417 /* make sure row indices within range */
1418 if (row < 0 || row >= n_row)
1420 stats [COLAMD_STATUS] = COLAMD_ERROR_row_index_out_of_bounds ;
1421 stats [COLAMD_INFO1] = col ;
1422 stats [COLAMD_INFO2] = row ;
1423 stats [COLAMD_INFO3] = n_row ;
1424 DEBUG0 ("colamd: row %d col %d out of bounds\n", row, col) ;
1425 return (FALSE) ;
1428 if (row <= last_row || Row [row].mark == col)
1430 /* row index are unsorted or repeated (or both), thus col */
1431 /* is jumbled. This is a notice, not an error condition. */
1432 stats [COLAMD_STATUS] = COLAMD_OK_BUT_JUMBLED ;
1433 stats [COLAMD_INFO1] = col ;
1434 stats [COLAMD_INFO2] = row ;
1435 (stats [COLAMD_INFO3]) ++ ;
1436 DEBUG1 ("colamd: row %d col %d unsorted/duplicate\n",row,col);
1439 if (Row [row].mark != col)
1441 Row [row].length++ ;
1443 else
1445 /* this is a repeated entry in the column, */
1446 /* it will be removed */
1447 Col [col].length-- ;
1450 /* mark the row as having been seen in this column */
1451 Row [row].mark = col ;
1453 last_row = row ;
1457 /* === Compute row pointers ========================================= */
1459 /* row form of the matrix starts directly after the column */
1460 /* form of matrix in A */
1461 Row [0].start = p [n_col] ;
1462 Row [0].p = Row [0].start ;
1463 Row [0].mark = -1 ;
1464 for (row = 1 ; row < n_row ; row++)
1466 Row [row].start = Row [row-1].start + Row [row-1].length ;
1467 Row [row].p = Row [row].start ;
1468 Row [row].mark = -1 ;
1471 /* === Create row form ============================================== */
1473 if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1475 /* if cols jumbled, watch for repeated row indices */
1476 for (col = 0 ; col < n_col ; col++)
1478 cp = A [p [col]] ;
1479 cp_end = A [p [col+1]] ;
1480 while (cp < cp_end)
1482 row = cp++ ;
1483 if (Row [row].mark != col)
1485 A [(Row [row].p)++] = col ;
1486 Row [row].mark = col ;
1491 else
1493 /* if cols not jumbled, we don't need the mark (this is faster) */
1494 for (col = 0 ; col < n_col ; col++)
1496 cp = A [p [col]] ;
1497 cp_end = A [p [col+1]] ;
1498 while (cp < cp_end)
1500 A [(Row [cp++].p)++] = col ;
1505 /* === Clear the row marks and set row degrees ====================== */
1507 for (row = 0 ; row < n_row ; row++)
1509 Row [row].mark = 0 ;
1510 Row [row].degree = Row [row].length ;
1513 /* === See if we need to re-create columns ========================== */
1515 if (stats [COLAMD_STATUS] == COLAMD_OK_BUT_JUMBLED)
1517 DEBUG0 ("colamd: reconstructing column form, matrix jumbled\n") ;
1519 if (!NDEBUG)
1521 /* make sure column lengths are correct */
1522 for (col = 0 ; col < n_col ; col++)
1524 p [col] = Col [col].length ;
1526 for (row = 0 ; row < n_row ; row++)
1528 rp = A [Row [row].start] ;
1529 rp_end = rp + Row [row].length ;
1530 while (rp < rp_end)
1532 p [rp++]-- ;
1535 for (col = 0 ; col < n_col ; col++)
1537 ASSERT (p [col] == 0) ;
1539 /* now p is all zero (different than when debugging is turned off) */
1540 } /* NDEBUG */
1542 /* === Compute col pointers ========================================= */
1544 /* col form of the matrix starts at A [0]. */
1545 /* Note, we may have a gap between the col form and the row */
1546 /* form if there were duplicate entries, if so, it will be */
1547 /* removed upon the first garbage collection */
1548 Col [0].start = 0 ;
1549 p [0] = Col [0].start ;
1550 for (col = 1 ; col < n_col ; col++)
1552 /* note that the lengths here are for pruned columns, i.e. */
1553 /* no duplicate row indices will exist for these columns */
1554 Col [col].start = Col [col-1].start + Col [col-1].length ;
1555 p [col] = Col [col].start ;
1558 /* === Re-create col form =========================================== */
1560 for (row = 0 ; row < n_row ; row++)
1562 rp = A [Row [row].start] ;
1563 rp_end = rp + Row [row].length ;
1564 while (rp < rp_end)
1566 A [(p [rp++])++] = row ;
1571 /* === Done. Matrix is not (or no longer) jumbled ================== */
1573 return (TRUE) ;
1577 /* ========================================================================== */
1578 /* === init_scoring ========================================================= */
1579 /* ========================================================================== */
1582 * Kills dense or empty columns and rows, calculates an initial score for
1583 * each column, and places all columns in the degree lists.init_rows_cols
1585 * @param n_row number of rows of A
1586 * @param n_col number of columns of A
1587 * @param Row of size n_row+1
1588 * @param Col of size n_col+1
1589 * @param A column form and row form of A
1590 * @param head of size n_col+1
1591 * @param knobs parameters
1592 * @param p_n_row2 size 1, number of non-dense, non-empty rows
1593 * @param p_n_col2 size 1, number of non-dense, non-empty columns
1594 * @param p_max_deg size 1, maximum row degree
1596 private static void init_scoring (int n_row, int n_col, Colamd_Row[] Row,
1597 Colamd_Col[] Col, int[] A, int[] head, double[] knobs,
1598 int[] p_n_row2, int[] p_n_col2, int[] p_max_deg)
1600 /* === Local variables ============================================== */
1602 int c ; /* a column index */
1603 int r, row ; /* a row index */
1604 int cp ; /* a column pointer */
1605 int deg ; /* degree of a row or column */
1606 int cp_end ; /* a pointer to the end of a column */
1607 int new_cp ; /* new column pointer */
1608 int col_length ; /* length of pruned column */
1609 int score ; /* current column score */
1610 int n_col2 ; /* number of non-dense, non-empty columns */
1611 int n_row2 ; /* number of non-dense, non-empty rows */
1612 int dense_row_count ; /* remove rows with more entries than this */
1613 int dense_col_count ; /* remove cols with more entries than this */
1614 int min_score ; /* smallest column score */
1615 int max_deg ; /* maximum row degree */
1616 int next_col ; /* Used to add to degree list.*/
1618 int debug_count = 0 ; /* debug only. */
1620 /* === Extract knobs ================================================ */
1622 /* Note: if knobs contains a NaN, this is undefined: */
1623 if (knobs [COLAMD_DENSE_ROW] < 0)
1625 /* only remove completely dense rows */
1626 dense_row_count = n_col-1 ;
1628 else
1630 dense_row_count = DENSE_DEGREE (knobs [COLAMD_DENSE_ROW], n_col) ;
1632 if (knobs [COLAMD_DENSE_COL] < 0)
1634 /* only remove completely dense columns */
1635 dense_col_count = n_row-1 ;
1637 else
1639 dense_col_count = DENSE_DEGREE (knobs [COLAMD_DENSE_COL],
1640 MIN (n_row, n_col)) ;
1643 DEBUG1 ("colamd: densecount: %d %d\n", dense_row_count, dense_col_count) ;
1644 max_deg = 0 ;
1645 n_col2 = n_col ;
1646 n_row2 = n_row ;
1648 /* === Kill empty columns =========================================== */
1650 /* Put the empty columns at the end in their natural order, so that LU */
1651 /* factorization can proceed as far as possible. */
1652 for (c = n_col-1 ; c >= 0 ; c--)
1654 deg = Col [c].length ;
1655 if (deg == 0)
1657 /* this is a empty column, kill and order it last */
1658 Col [c].order = --n_col2 ;
1659 KILL_PRINCIPAL_COL (Col, c) ;
1662 DEBUG1 ("colamd: null columns killed: %d\n", n_col - n_col2) ;
1664 /* === Kill dense columns =========================================== */
1666 /* Put the dense columns at the end, in their natural order */
1667 for (c = n_col-1 ; c >= 0 ; c--)
1669 /* skip any dead columns */
1670 if (COL_IS_DEAD (Col, c))
1672 continue ;
1674 deg = Col [c].length ;
1675 if (deg > dense_col_count)
1677 /* this is a dense column, kill and order it last */
1678 Col [c].order = --n_col2 ;
1679 /* decrement the row degrees */
1680 cp = A [Col [c].start] ;
1681 cp_end = cp + Col [c].length ;
1682 while (cp < cp_end)
1684 Row [cp++].degree-- ;
1686 KILL_PRINCIPAL_COL (Col, c) ;
1689 DEBUG1 ("colamd: Dense and null columns killed: %d\n", n_col - n_col2) ;
1691 /* === Kill dense and empty rows ==================================== */
1693 for (r = 0 ; r < n_row ; r++)
1695 deg = Row [r].degree ;
1696 ASSERT (deg >= 0 && deg <= n_col) ;
1697 if (deg > dense_row_count || deg == 0)
1699 /* kill a dense or empty row */
1700 KILL_ROW (Row, r) ;
1701 --n_row2 ;
1703 else
1705 /* keep track of max degree of remaining rows */
1706 max_deg = MAX (max_deg, deg) ;
1709 DEBUG1 ("colamd: Dense and null rows killed: %d\n", n_row - n_row2) ;
1711 /* === Compute initial column scores ================================ */
1713 /* At this point the row degrees are accurate. They reflect the number */
1714 /* of "live" (non-dense) columns in each row. No empty rows exist. */
1715 /* Some "live" columns may contain only dead rows, however. These are */
1716 /* pruned in the code below. */
1718 /* now find the initial matlab score for each column */
1719 for (c = n_col-1 ; c >= 0 ; c--)
1721 /* skip dead column */
1722 if (COL_IS_DEAD (Col, c))
1724 continue ;
1726 score = 0 ;
1727 cp = A [Col [c].start] ;
1728 new_cp = cp ;
1729 cp_end = cp + Col [c].length ;
1730 while (cp < cp_end)
1732 /* get a row */
1733 row = cp++ ;
1734 /* skip if dead */
1735 if (ROW_IS_DEAD (Row, row))
1737 continue ;
1739 /* FIXME: compact the column */
1740 new_cp = row ;
1741 // new_cp++ = row ;
1742 /* add row's external degree */
1743 score += Row [row].degree - 1 ;
1744 /* guard against integer overflow */
1745 score = MIN (score, n_col) ;
1747 /* determine pruned column length */
1748 col_length = (new_cp - A [Col [c].start]) ;
1749 if (col_length == 0)
1751 /* a newly-made null column (all rows in this col are "dense" */
1752 /* and have already been killed) */
1753 DEBUG2 ("Newly null killed: %d\n", c) ;
1754 Col [c].order = --n_col2 ;
1755 KILL_PRINCIPAL_COL (Col, c) ;
1757 else
1759 /* set column length and set score */
1760 ASSERT (score >= 0) ;
1761 ASSERT (score <= n_col) ;
1762 Col [c].length = col_length ;
1763 Col [c].score = score ;
1766 DEBUG1 ("colamd: Dense, null, and newly-null columns killed: %d\n",
1767 n_col-n_col2) ;
1769 /* At this point, all empty rows and columns are dead. All live columns */
1770 /* are "clean" (containing no dead rows) and simplicial (no supercolumns */
1771 /* yet). Rows may contain dead columns, but all live rows contain at */
1772 /* least one live column. */
1774 if (!NDEBUG)
1776 debug_structures (n_row, n_col, Row, Col, A, n_col2) ;
1779 /* === Initialize degree lists ========================================== */
1781 if (!NDEBUG)
1783 debug_count = 0 ;
1786 /* clear the hash buckets */
1787 for (c = 0 ; c <= n_col ; c++)
1789 head [c] = EMPTY ;
1791 min_score = n_col ;
1792 /* place in reverse order, so low column indices are at the front */
1793 /* of the lists. This is to encourage natural tie-breaking */
1794 for (c = n_col-1 ; c >= 0 ; c--)
1796 /* only add principal columns to degree lists */
1797 if (COL_IS_ALIVE (Col, c))
1799 DEBUG4 ("place %d score %d minscore %d ncol %d\n",
1800 c, Col [c].score, min_score, n_col) ;
1802 /* === Add columns score to DList =============================== */
1804 score = Col [c].score ;
1806 ASSERT (min_score >= 0) ;
1807 ASSERT (min_score <= n_col) ;
1808 ASSERT (score >= 0) ;
1809 ASSERT (score <= n_col) ;
1810 ASSERT (head [score] >= EMPTY) ;
1812 /* now add this column to dList at proper score location */
1813 next_col = head [score] ;
1814 Col [c].prev = EMPTY ;
1815 Col [c].degree_next = next_col ;
1817 /* if there already was a column with the same score, set its */
1818 /* previous pointer to this new column */
1819 if (next_col != EMPTY)
1821 Col [next_col].prev = c ;
1823 head [score] = c ;
1825 /* see if this score is less than current min */
1826 min_score = MIN (min_score, score) ;
1828 if (!NDEBUG)
1830 debug_count++ ;
1836 if (!NDEBUG)
1838 DEBUG1 ("colamd: Live cols %d out of %d, non-princ: %d\n",
1839 debug_count, n_col, n_col-debug_count) ;
1840 ASSERT (debug_count == n_col2) ;
1841 debug_deg_lists (n_row, n_col, Row, Col, head, min_score, n_col2, max_deg) ;
1842 } /* NDEBUG */
1844 /* === Return number of remaining columns, and max row degree ======= */
1846 p_n_col2[0] = n_col2 ;
1847 p_n_row2[0] = n_row2 ;
1848 p_max_deg[0] = max_deg ;
1852 /* ========================================================================== */
1853 /* === find_ordering ======================================================== */
1854 /* ========================================================================== */
1857 * Order the principal columns of the supercolumn form of the matrix
1858 * (no supercolumns on input). Uses a minimum approximate column minimum
1859 * degree ordering method. Not user-callable.
1861 * @param n_row number of rows of A
1862 * @param n_col number of columns of A
1863 * @param Alen size of A, 2*nnz + n_col or larger
1864 * @param Row of size n_row+1
1865 * @param Col of size n_col+1
1866 * @param A column form and row form of A
1867 * @param head of size n_col+1
1868 * @param n_col2 Remaining columns to order
1869 * @param max_deg Maximum row degree
1870 * @param pfree index of first free slot (2*nnz on entry)
1871 * @param aggressive
1872 * @return the number of garbage collections
1874 private static int find_ordering (int n_row, int n_col, int Alen,
1875 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int[] head, int n_col2,
1876 int max_deg, int pfree, int aggressive)
1878 /* === Local variables ============================================== */
1880 int k ; /* current pivot ordering step */
1881 int pivot_col ; /* current pivot column */
1882 int cp ; /* a column pointer */
1883 int rp ; /* a row pointer */
1884 int pivot_row ; /* current pivot row */
1885 int new_cp ; /* modified column pointer */
1886 int new_rp ; /* modified row pointer */
1887 int pivot_row_start ; /* pointer to start of pivot row */
1888 int pivot_row_degree ; /* number of columns in pivot row */
1889 int pivot_row_length ; /* number of supercolumns in pivot row */
1890 int pivot_col_score ; /* score of pivot column */
1891 int needed_memory ; /* free space needed for pivot row */
1892 int cp_end ; /* pointer to the end of a column */
1893 int rp_end ; /* pointer to the end of a row */
1894 int row ; /* a row index */
1895 int col ; /* a column index */
1896 int max_score ; /* maximum possible score */
1897 int cur_score ; /* score of current column */
1898 /*FIXME: unsigned*/ int hash ; /* hash value for supernode detection */
1899 int head_column ; /* head of hash bucket */
1900 int first_col ; /* first column in hash bucket */
1901 int tag_mark ; /* marker value for mark array */
1902 int row_mark ; /* Row [row].shared2.mark */
1903 int set_difference ; /* set difference size of row with pivot row */
1904 int min_score ; /* smallest column score */
1905 int col_thickness ; /* "thickness" (no. of columns in a supercol) */
1906 int max_mark ; /* maximum value of tag_mark */
1907 int pivot_col_thickness ; /* number of columns represented by pivot col */
1908 int prev_col ; /* Used by Dlist operations. */
1909 int next_col ; /* Used by Dlist operations. */
1910 int ngarbage ; /* number of garbage collections performed */
1912 int debug_d ; /* debug loop counter */
1913 int debug_step = 0 ; /* debug loop counter */
1915 /* === Initialization and clear mark ================================ */
1917 max_mark = Int_MAX - n_col ;
1918 tag_mark = clear_mark (0, max_mark, n_row, Row) ;
1919 min_score = 0 ;
1920 ngarbage = 0 ;
1921 DEBUG1 ("colamd: Ordering, n_col2=%d\n", n_col2) ;
1923 /* === Order the columns ============================================ */
1925 for (k = 0 ; k < n_col2 ; /* 'k' is incremented below */)
1928 if (!NDEBUG)
1930 if (debug_step % 100 == 0)
1932 DEBUG2 ("\n... Step k: %d out of n_col2: %d\n", k, n_col2) ;
1934 else
1936 DEBUG3 ("\n----------Step k: %d out of n_col2: %d\n", k, n_col2) ;
1938 debug_step++ ;
1939 debug_deg_lists (n_row, n_col, Row, Col, head,
1940 min_score, n_col2-k, max_deg) ;
1941 debug_matrix (n_row, n_col, Row, Col, A) ;
1942 } /* NDEBUG */
1944 /* === Select pivot column, and order it ============================ */
1946 /* make sure degree list isn't empty */
1947 ASSERT (min_score >= 0) ;
1948 ASSERT (min_score <= n_col) ;
1949 ASSERT (head [min_score] >= EMPTY) ;
1951 if (!NDEBUG)
1953 for (debug_d = 0 ; debug_d < min_score ; debug_d++)
1955 ASSERT (head [debug_d] == EMPTY) ;
1957 } /* NDEBUG */
1959 /* get pivot column from head of minimum degree list */
1960 while (head [min_score] == EMPTY && min_score < n_col)
1962 min_score++ ;
1964 pivot_col = head [min_score] ;
1965 ASSERT (pivot_col >= 0 && pivot_col <= n_col) ;
1966 next_col = Col [pivot_col].degree_next ;
1967 head [min_score] = next_col ;
1968 if (next_col != EMPTY)
1970 Col [next_col].prev = EMPTY ;
1973 ASSERT (COL_IS_ALIVE (Col, pivot_col)) ;
1975 /* remember score for defrag check */
1976 pivot_col_score = Col [pivot_col].score ;
1978 /* the pivot column is the kth column in the pivot order */
1979 Col [pivot_col].order = k ;
1981 /* increment order count by column thickness */
1982 pivot_col_thickness = Col [pivot_col].thickness ;
1983 k += pivot_col_thickness ;
1984 ASSERT (pivot_col_thickness > 0) ;
1985 DEBUG3 ("Pivot col: %d thick %d\n", pivot_col, pivot_col_thickness) ;
1987 /* === Garbage_collection, if necessary ============================= */
1989 needed_memory = MIN (pivot_col_score, n_col - k) ;
1990 if (pfree + needed_memory >= Alen)
1992 pfree = garbage_collection (n_row, n_col, Row, Col, A, A [pfree]) ;
1993 ngarbage++ ;
1994 /* after garbage collection we will have enough */
1995 ASSERT (pfree + needed_memory < Alen) ;
1996 /* garbage collection has wiped out the Row[].shared2.mark array */
1997 tag_mark = clear_mark (0, max_mark, n_row, Row) ;
1999 if (!NDEBUG)
2001 debug_matrix (n_row, n_col, Row, Col, A) ;
2002 } /* NDEBUG */
2005 /* === Compute pivot row pattern ==================================== */
2007 /* get starting location for this new merged row */
2008 pivot_row_start = pfree ;
2010 /* initialize new row counts to zero */
2011 pivot_row_degree = 0 ;
2013 /* tag pivot column as having been visited so it isn't included */
2014 /* in merged pivot row */
2015 Col [pivot_col].thickness = -pivot_col_thickness ;
2017 /* pivot row is the union of all rows in the pivot column pattern */
2018 cp = A [Col [pivot_col].start] ;
2019 cp_end = cp + Col [pivot_col].length ;
2020 while (cp < cp_end)
2022 /* get a row */
2023 row = cp++ ;
2024 DEBUG4 ("Pivot col pattern %d %d\n", ROW_IS_ALIVE (Row, row), row) ;
2025 /* skip if row is dead */
2026 if (ROW_IS_ALIVE (Row, row))
2028 rp = A [Row [row].start] ;
2029 rp_end = rp + Row [row].length ;
2030 while (rp < rp_end)
2032 /* get a column */
2033 col = rp++ ;
2034 /* add the column, if alive and untagged */
2035 col_thickness = Col [col].thickness ;
2036 if (col_thickness > 0 && COL_IS_ALIVE (Col, col))
2038 /* tag column in pivot row */
2039 Col [col].thickness = -col_thickness ;
2040 ASSERT (pfree < Alen) ;
2041 /* place column in pivot row */
2042 A [pfree++] = col ;
2043 pivot_row_degree += col_thickness ;
2049 /* clear tag on pivot column */
2050 Col [pivot_col].thickness = pivot_col_thickness ;
2051 max_deg = MAX (max_deg, pivot_row_degree) ;
2053 if (!NDEBUG)
2055 DEBUG3 ("check2\n") ;
2056 debug_mark (n_row, Row, tag_mark, max_mark) ;
2057 } /* NDEBUG */
2059 /* === Kill all rows used to construct pivot row ==================== */
2061 /* also kill pivot row, temporarily */
2062 cp = A [Col [pivot_col].start] ;
2063 cp_end = cp + Col [pivot_col].length ;
2064 while (cp < cp_end)
2066 /* may be killing an already dead row */
2067 row = cp++ ;
2068 DEBUG3 ("Kill row in pivot col: %d\n", row) ;
2069 KILL_ROW (Row, row) ;
2072 /* === Select a row index to use as the new pivot row =============== */
2074 pivot_row_length = pfree - pivot_row_start ;
2075 if (pivot_row_length > 0)
2077 /* pick the "pivot" row arbitrarily (first row in col) */
2078 pivot_row = A [Col [pivot_col].start] ;
2079 DEBUG3 ("Pivotal row is %d\n", pivot_row) ;
2081 else
2083 /* there is no pivot row, since it is of zero length */
2084 pivot_row = EMPTY ;
2085 ASSERT (pivot_row_length == 0) ;
2087 ASSERT (Col [pivot_col].length > 0 || pivot_row_length == 0) ;
2089 /* === Approximate degree computation =============================== */
2091 /* Here begins the computation of the approximate degree. The column */
2092 /* score is the sum of the pivot row "length", plus the size of the */
2093 /* set differences of each row in the column minus the pattern of the */
2094 /* pivot row itself. The column ("thickness") itself is also */
2095 /* excluded from the column score (we thus use an approximate */
2096 /* external degree). */
2098 /* The time taken by the following code (compute set differences, and */
2099 /* add them up) is proportional to the size of the data structure */
2100 /* being scanned - that is, the sum of the sizes of each column in */
2101 /* the pivot row. Thus, the amortized time to compute a column score */
2102 /* is proportional to the size of that column (where size, in this */
2103 /* context, is the column "length", or the number of row indices */
2104 /* in that column). The number of row indices in a column is */
2105 /* monotonically non-decreasing, from the length of the original */
2106 /* column on input to colamd. */
2108 /* === Compute set differences ====================================== */
2110 DEBUG3 ("** Computing set differences phase. **\n") ;
2112 /* pivot row is currently dead - it will be revived later. */
2114 DEBUG3 ("Pivot row: ") ;
2115 /* for each column in pivot row */
2116 rp = A [pivot_row_start] ;
2117 rp_end = rp + pivot_row_length ;
2118 while (rp < rp_end)
2120 col = rp++ ;
2121 ASSERT (COL_IS_ALIVE (Col, col) && col != pivot_col) ;
2122 DEBUG3 ("Col: %d\n", col) ;
2124 /* clear tags used to construct pivot row pattern */
2125 col_thickness = -Col [col].thickness ;
2126 ASSERT (col_thickness > 0) ;
2127 Col [col].thickness = col_thickness ;
2129 /* === Remove column from degree list =========================== */
2131 cur_score = Col [col].score ;
2132 prev_col = Col [col].prev ;
2133 next_col = Col [col].degree_next ;
2134 ASSERT (cur_score >= 0) ;
2135 ASSERT (cur_score <= n_col) ;
2136 ASSERT (cur_score >= EMPTY) ;
2137 if (prev_col == EMPTY)
2139 head [cur_score] = next_col ;
2141 else
2143 Col [prev_col].degree_next = next_col ;
2145 if (next_col != EMPTY)
2147 Col [next_col].prev = prev_col ;
2150 /* === Scan the column ========================================== */
2152 cp = A [Col [col].start] ;
2153 cp_end = cp + Col [col].length ;
2154 while (cp < cp_end)
2156 /* get a row */
2157 row = cp++ ;
2158 row_mark = Row [row].mark ;
2159 /* skip if dead */
2160 if (ROW_IS_MARKED_DEAD (row_mark))
2162 continue ;
2164 ASSERT (row != pivot_row) ;
2165 set_difference = row_mark - tag_mark ;
2166 /* check if the row has been seen yet */
2167 if (set_difference < 0)
2169 ASSERT (Row [row].degree <= max_deg) ;
2170 set_difference = Row [row].degree ;
2172 /* subtract column thickness from this row's set difference */
2173 set_difference -= col_thickness ;
2174 ASSERT (set_difference >= 0) ;
2175 /* absorb this row if the set difference becomes zero */
2176 if (set_difference == 0 && aggressive != 0)
2178 DEBUG3 ("aggressive absorption. Row: %d\n", row) ;
2179 KILL_ROW (Row, row) ;
2181 else
2183 /* save the new mark */
2184 Row [row].mark = set_difference + tag_mark ;
2189 if (!NDEBUG)
2191 debug_deg_lists (n_row, n_col, Row, Col, head,
2192 min_score, n_col2-k-pivot_row_degree, max_deg) ;
2193 } /* NDEBUG */
2195 /* === Add up set differences for each column ======================= */
2197 DEBUG3 ("** Adding set differences phase. **\n") ;
2199 /* for each column in pivot row */
2200 rp = A [pivot_row_start] ;
2201 rp_end = rp + pivot_row_length ;
2202 while (rp < rp_end)
2204 /* get a column */
2205 col = rp++ ;
2206 ASSERT (COL_IS_ALIVE (Col, col) && col != pivot_col) ;
2207 hash = 0 ;
2208 cur_score = 0 ;
2209 cp = A [Col [col].start] ;
2210 /* compact the column */
2211 new_cp = cp ;
2212 cp_end = cp + Col [col].length ;
2214 DEBUG4 ("Adding set diffs for Col: %d.\n", col) ;
2216 while (cp < cp_end)
2218 /* get a row */
2219 row = cp++ ;
2220 ASSERT(row >= 0 && row < n_row) ;
2221 row_mark = Row [row].mark ;
2222 /* skip if dead */
2223 if (ROW_IS_MARKED_DEAD (row_mark))
2225 DEBUG4 (" Row %d, dead\n", row) ;
2226 continue ;
2228 DEBUG4 (" Row %d, set diff %d\n", row, row_mark-tag_mark);
2229 ASSERT (row_mark >= tag_mark) ;
2230 /* FIXME: compact the column */
2231 new_cp = row ;
2232 // new_cp++ = row ;
2233 /* compute hash function */
2234 hash += row ;
2235 /* add set difference */
2236 cur_score += row_mark - tag_mark ;
2237 /* integer overflow... */
2238 cur_score = MIN (cur_score, n_col) ;
2241 /* recompute the column's length */
2242 Col [col].length = (int) (new_cp - A [Col [col].start]) ;
2244 /* === Further mass elimination ================================= */
2246 if (Col [col].length == 0)
2248 DEBUG4 ("further mass elimination. Col: %d\n", col) ;
2249 /* nothing left but the pivot row in this column */
2250 KILL_PRINCIPAL_COL (Col, col) ;
2251 pivot_row_degree -= Col [col].thickness ;
2252 ASSERT (pivot_row_degree >= 0) ;
2253 /* order it */
2254 Col [col].order = k ;
2255 /* increment order count by column thickness */
2256 k += Col [col].thickness ;
2258 else
2260 /* === Prepare for supercolumn detection ==================== */
2262 DEBUG4 ("Preparing supercol detection for Col: %d.\n", col) ;
2264 /* save score so far */
2265 Col [col].score = cur_score ;
2267 /* add column to hash table, for supercolumn detection */
2268 hash %= n_col + 1 ;
2270 DEBUG4 (" Hash = %d, n_col = %d.\n", hash, n_col) ;
2271 ASSERT (((int) hash) <= n_col) ;
2273 head_column = head [hash] ;
2274 if (head_column > EMPTY)
2276 /* degree list "hash" is non-empty, use prev (shared3) of */
2277 /* first column in degree list as head of hash bucket */
2278 first_col = Col [head_column].headhash ;
2279 Col [head_column].headhash = col ;
2281 else
2283 /* degree list "hash" is empty, use head as hash bucket */
2284 first_col = - (head_column + 2) ;
2285 head [hash] = - (col + 2) ;
2287 Col [col].hash_next = first_col ;
2289 /* save hash function in Col [col].shared3.hash */
2290 Col [col].hash = (int) hash ;
2291 ASSERT (COL_IS_ALIVE (Col, col)) ;
2295 /* The approximate external column degree is now computed. */
2297 /* === Supercolumn detection ======================================== */
2299 DEBUG3 ("** Supercolumn detection phase. **\n") ;
2301 if (!NDEBUG)
2303 detect_super_cols (n_col, Row,
2304 Col, A, head, pivot_row_start, pivot_row_length) ;
2305 } else {
2306 detect_super_cols (Col, A, head, pivot_row_start, pivot_row_length) ;
2309 /* === Kill the pivotal column ====================================== */
2311 KILL_PRINCIPAL_COL (Col, pivot_col) ;
2313 /* === Clear mark =================================================== */
2315 tag_mark = clear_mark (tag_mark+max_deg+1, max_mark, n_row, Row) ;
2317 if (!NDEBUG)
2319 DEBUG3 ("check3\n") ;
2320 debug_mark (n_row, Row, tag_mark, max_mark) ;
2321 } /* NDEBUG */
2323 /* === Finalize the new pivot row, and column scores ================ */
2325 DEBUG3 ("** Finalize scores phase. **\n") ;
2327 /* for each column in pivot row */
2328 rp = A [pivot_row_start] ;
2329 /* compact the pivot row */
2330 new_rp = rp ;
2331 rp_end = rp + pivot_row_length ;
2332 while (rp < rp_end)
2334 col = rp++ ;
2335 /* skip dead columns */
2336 if (COL_IS_DEAD (Col, col))
2338 continue ;
2340 new_rp = col ; // FIXME: pointer
2341 //new_rp++ = col ;
2342 /* add new pivot row to column */
2343 A [Col [col].start + (Col [col].length++)] = pivot_row ;
2345 /* retrieve score so far and add on pivot row's degree. */
2346 /* (we wait until here for this in case the pivot */
2347 /* row's degree was reduced due to mass elimination). */
2348 cur_score = Col [col].score + pivot_row_degree ;
2350 /* calculate the max possible score as the number of */
2351 /* external columns minus the 'k' value minus the */
2352 /* columns thickness */
2353 max_score = n_col - k - Col [col].thickness ;
2355 /* make the score the external degree of the union-of-rows */
2356 cur_score -= Col [col].thickness ;
2358 /* make sure score is less or equal than the max score */
2359 cur_score = MIN (cur_score, max_score) ;
2360 ASSERT (cur_score >= 0) ;
2362 /* store updated score */
2363 Col [col].score = cur_score ;
2365 /* === Place column back in degree list ========================= */
2367 ASSERT (min_score >= 0) ;
2368 ASSERT (min_score <= n_col) ;
2369 ASSERT (cur_score >= 0) ;
2370 ASSERT (cur_score <= n_col) ;
2371 ASSERT (head [cur_score] >= EMPTY) ;
2372 next_col = head [cur_score] ;
2373 Col [col].degree_next = next_col ;
2374 Col [col].prev = EMPTY ;
2375 if (next_col != EMPTY)
2377 Col [next_col].prev = col ;
2379 head [cur_score] = col ;
2381 /* see if this score is less than current min */
2382 min_score = MIN (min_score, cur_score) ;
2386 if (!NDEBUG)
2388 debug_deg_lists (n_row, n_col, Row, Col, head,
2389 min_score, n_col2-k, max_deg) ;
2390 } /* NDEBUG */
2392 /* === Resurrect the new pivot row ================================== */
2394 if (pivot_row_degree > 0)
2396 /* update pivot row length to reflect any cols that were killed */
2397 /* during super-col detection and mass elimination */
2398 Row [pivot_row].start = pivot_row_start ;
2399 Row [pivot_row].length = (int) (new_rp - A [pivot_row_start]) ;
2400 ASSERT (Row [pivot_row].length > 0) ;
2401 Row [pivot_row].degree = pivot_row_degree ;
2402 Row [pivot_row].mark = 0 ;
2403 /* pivot row is no longer dead */
2405 DEBUG1 ("Resurrect Pivot_row %d deg: %d\n",
2406 pivot_row, pivot_row_degree) ;
2410 /* === All principal columns have now been ordered ================== */
2412 return (ngarbage) ;
2416 /* ========================================================================== */
2417 /* === order_children ======================================================= */
2418 /* ========================================================================== */
2421 * The find_ordering routine has ordered all of the principal columns (the
2422 * representatives of the supercolumns). The non-principal columns have not
2423 * yet been ordered. This routine orders those columns by walking up the
2424 * parent tree (a column is a child of the column which absorbed it). The
2425 * final permutation vector is then placed in p [0 ... n_col-1], with p [0]
2426 * being the first column, and p [n_col-1] being the last. It doesn't look
2427 * like it at first glance, but be assured that this routine takes time linear
2428 * in the number of columns. Although not immediately obvious, the time
2429 * taken by this routine is O (n_col), that is, linear in the number of
2430 * columns. Not user-callable.
2432 * @param n_col number of columns of A
2433 * @param Col of size n_col+1
2434 * @param p p [0 ... n_col-1] is the column permutation
2436 private static void order_children (int n_col, Colamd_Col[] Col, int[] p)
2438 /* === Local variables ============================================== */
2440 int i ; /* loop counter for all columns */
2441 int c ; /* column index */
2442 int parent ; /* index of column's parent */
2443 int order ; /* column's order */
2445 /* === Order each non-principal column ============================== */
2447 for (i = 0 ; i < n_col ; i++)
2449 /* find an un-ordered non-principal column */
2450 ASSERT (COL_IS_DEAD (Col, i)) ;
2451 if (!COL_IS_DEAD_PRINCIPAL (Col, i) && Col [i].order == EMPTY)
2453 parent = i ;
2454 /* once found, find its principal parent */
2457 parent = Col [parent].parent ;
2458 } while (!COL_IS_DEAD_PRINCIPAL (Col, parent)) ;
2460 /* now, order all un-ordered non-principal columns along path */
2461 /* to this parent. collapse tree at the same time */
2462 c = i ;
2463 /* get order of parent */
2464 order = Col [parent].order ;
2468 ASSERT (Col [c].order == EMPTY) ;
2470 /* order this column */
2471 Col [c].order = order++ ;
2472 /* collaps tree */
2473 Col [c].parent = parent ;
2475 /* get immediate parent of this column */
2476 c = Col [c].parent ;
2478 /* continue until we hit an ordered column. There are */
2479 /* guarranteed not to be anymore unordered columns */
2480 /* above an ordered column */
2481 } while (Col [c].order == EMPTY) ;
2483 /* re-order the super_col parent to largest order for this group */
2484 Col [parent].order = order ;
2488 /* === Generate the permutation ===================================== */
2490 for (c = 0 ; c < n_col ; c++)
2492 p [Col [c].order] = c ;
2497 /* ========================================================================== */
2498 /* === detect_super_cols ==================================================== */
2499 /* ========================================================================== */
2502 * Detects supercolumns by finding matches between columns in the hash buckets.
2503 * Check amongst columns in the set A [row_start ... row_start + row_length-1].
2504 * The columns under consideration are currently *not* in the degree lists,
2505 * and have already been placed in the hash buckets.
2507 * The hash bucket for columns whose hash function is equal to h is stored
2508 * as follows:
2510 * if head [h] is >= 0, then head [h] contains a degree list, so:
2512 * head [h] is the first column in degree bucket h.
2513 * Col [head [h]].headhash gives the first column in hash bucket h.
2515 * otherwise, the degree list is empty, and:
2517 * -(head [h] + 2) is the first column in hash bucket h.
2519 * For a column c in a hash bucket, Col [c].shared3.prev is NOT a "previous
2520 * column" pointer. Col [c].shared3.hash is used instead as the hash number
2521 * for that column. The value of Col [c].shared4.hash_next is the next column
2522 * in the same hash bucket.
2524 * Assuming no, or "few" hash collisions, the time taken by this routine is
2525 * linear in the sum of the sizes (lengths) of each column whose score has
2526 * just been computed in the approximate degree computation.
2527 * Not user-callable.
2529 * @param Col of size n_col+1
2530 * @param A row indices of A
2531 * @param head head of degree lists and hash buckets
2532 * @param row_start pointer to set of columns to check
2533 * @param row_length number of columns to check
2535 private static void detect_super_cols (Colamd_Col[] Col, int[] A,
2536 int[] head, int row_start, int row_length)
2538 detect_super_cols (0, null, Col, A, head, row_start, row_length) ;
2542 * debug version
2544 * @param n_col number of columns of A
2545 * @param Row of size n_row+1
2546 * @param Col of size n_col+1
2547 * @param A row indices of A
2548 * @param head head of degree lists and hash buckets
2549 * @param row_start pointer to set of columns to check
2550 * @param row_length number of columns to check
2552 private static void detect_super_cols (int n_col, Colamd_Row[] Row,
2553 Colamd_Col[] Col, int[] A, int[] head, int row_start, int row_length)
2555 /* === Local variables ============================================== */
2557 int hash ; /* hash value for a column */
2558 int rp ; /* pointer to a row */
2559 int c ; /* a column index */
2560 int super_c ; /* column index of the column to absorb into */
2561 int cp1 ; /* column pointer for column super_c */
2562 int cp2 ; /* column pointer for column c */
2563 int length ; /* length of column super_c */
2564 int prev_c ; /* column preceding c in hash bucket */
2565 int i ; /* loop counter */
2566 int rp_end ; /* pointer to the end of the row */
2567 int col ; /* a column index in the row to check */
2568 int head_column ; /* first column in hash bucket or degree list */
2569 int first_col ; /* first column in hash bucket */
2571 /* === Consider each column in the row ============================== */
2573 rp = A [row_start] ;
2574 rp_end = rp + row_length ;
2575 while (rp < rp_end)
2577 col = rp++ ;
2578 if (COL_IS_DEAD (Col, col))
2580 continue ;
2583 /* get hash number for this column */
2584 hash = Col [col].hash ;
2585 ASSERT (hash <= n_col) ;
2587 /* === Get the first column in this hash bucket ===================== */
2589 head_column = head [hash] ;
2590 if (head_column > EMPTY)
2592 first_col = Col [head_column].headhash ;
2594 else
2596 first_col = - (head_column + 2) ;
2599 /* === Consider each column in the hash bucket ====================== */
2601 for (super_c = first_col ; super_c != EMPTY ;
2602 super_c = Col [super_c].hash_next)
2604 ASSERT (COL_IS_ALIVE (Col, super_c)) ;
2605 ASSERT (Col [super_c].hash == hash) ;
2606 length = Col [super_c].length ;
2608 /* prev_c is the column preceding column c in the hash bucket */
2609 prev_c = super_c ;
2611 /* === Compare super_c with all columns after it ================ */
2613 for (c = Col [super_c].hash_next ;
2614 c != EMPTY ; c = Col [c].hash_next)
2616 ASSERT (c != super_c) ;
2617 ASSERT (COL_IS_ALIVE (Col, c)) ;
2618 ASSERT (Col [c].hash == hash) ;
2620 /* not identical if lengths or scores are different */
2621 if (Col [c].length != length ||
2622 Col [c].score != Col [super_c].score)
2624 prev_c = c ;
2625 continue ;
2628 /* compare the two columns */
2629 cp1 = A [Col [super_c].start] ;
2630 cp2 = A [Col [c].start] ;
2632 for (i = 0 ; i < length ; i++)
2634 /* the columns are "clean" (no dead rows) */
2635 ASSERT (ROW_IS_ALIVE (Row, cp1)) ;
2636 ASSERT (ROW_IS_ALIVE (Row, cp2)) ;
2637 /* row indices will same order for both supercols, */
2638 /* no gather scatter nessasary */
2639 if (cp1++ != cp2++)
2641 break ;
2645 /* the two columns are different if the for-loop "broke" */
2646 if (i != length)
2648 prev_c = c ;
2649 continue ;
2652 /* === Got it! two columns are identical =================== */
2654 ASSERT (Col [c].score == Col [super_c].score) ;
2656 Col [super_c].thickness += Col [c].thickness ;
2657 Col [c].parent = super_c ;
2658 KILL_NON_PRINCIPAL_COL (Col, c) ;
2659 /* order c later, in order_children() */
2660 Col [c].order = EMPTY ;
2661 /* remove c from hash bucket */
2662 Col [prev_c].hash_next = Col [c].hash_next ;
2666 /* === Empty this hash bucket ======================================= */
2668 if (head_column > EMPTY)
2670 /* corresponding degree list "hash" is not empty */
2671 Col [head_column].headhash = EMPTY ;
2673 else
2675 /* corresponding degree list "hash" is empty */
2676 head [hash] = EMPTY ;
2682 /* ========================================================================== */
2683 /* === garbage_collection =================================================== */
2684 /* ========================================================================== */
2687 * Defragments and compacts columns and rows in the workspace A. Used when
2688 * all avaliable memory has been used while performing row merging. Returns
2689 * the index of the first free position in A, after garbage collection. The
2690 * time taken by this routine is linear in the size of the array A, which is
2691 * itself linear in the number of nonzeros in the input matrix.
2692 * Not user-callable.
2694 * @param n_row number of rows
2695 * @param n_col number of columns
2696 * @param Row row info
2697 * @param Col column info
2698 * @param A A [0 ... Alen-1] holds the matrix
2699 * @param pfree
2700 * @return the new value of pfree
2702 private static int garbage_collection (int n_row, int n_col,
2703 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int pfree)
2705 /* === Local variables ============================================== */
2707 int psrc ; /* source pointer */
2708 int pdest ; /* destination pointer */
2709 int j ; /* counter */
2710 int r ; /* a row index */
2711 int c ; /* a column index */
2712 int length ; /* length of a row or column */
2714 int debug_rows = 0 ;
2715 if (!NDEBUG)
2717 DEBUG2 ("Defrag..\n") ;
2718 for (psrc = A[0] ; psrc < pfree ; psrc++) ASSERT (psrc >= 0) ;
2719 debug_rows = 0 ;
2722 /* === Defragment the columns ======================================= */
2724 pdest = A[0] ;
2725 for (c = 0 ; c < n_col ; c++)
2727 if (COL_IS_ALIVE (Col, c))
2729 psrc = A [Col [c].start] ;
2731 /* move and compact the column */
2732 ASSERT (pdest <= psrc) ;
2733 Col [c].start = (int) (pdest - A [0]) ;
2734 length = Col [c].length ;
2735 for (j = 0 ; j < length ; j++)
2737 r = psrc++ ;
2738 if (ROW_IS_ALIVE (Row, r))
2740 pdest = r ; // FIXME: pointer
2741 //pdest++ = r ;
2744 Col [c].length = (int) (pdest - A [Col [c].start]) ;
2748 /* === Prepare to defragment the rows =============================== */
2750 for (r = 0 ; r < n_row ; r++)
2752 if (ROW_IS_DEAD (Row, r) || (Row [r].length == 0))
2754 /* This row is already dead, or is of zero length. Cannot compact
2755 * a row of zero length, so kill it. NOTE: in the current version,
2756 * there are no zero-length live rows. Kill the row (for the first
2757 * time, or again) just to be safe. */
2758 KILL_ROW (Row, r) ;
2760 else
2762 /* save first column index in Row [r].shared2.first_column */
2763 psrc = A [Row [r].start] ;
2764 Row [r].first_column = psrc ;
2765 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2766 /* flag the start of the row with the one's complement of row */
2767 psrc = ONES_COMPLEMENT (r) ;
2768 if (!NDEBUG)
2770 debug_rows++ ;
2771 } /* NDEBUG */
2775 /* === Defragment the rows ========================================== */
2777 psrc = pdest ;
2778 while (psrc < pfree)
2780 /* find a negative number ... the start of a row */
2781 if (psrc++ < 0)
2783 psrc-- ;
2784 /* get the row index */
2785 r = ONES_COMPLEMENT (psrc) ;
2786 ASSERT (r >= 0 && r < n_row) ;
2787 /* restore first column index */
2788 psrc = Row [r].first_column ;
2789 ASSERT (ROW_IS_ALIVE (Row, r)) ;
2790 ASSERT (Row [r].length > 0) ;
2791 /* move and compact the row */
2792 ASSERT (pdest <= psrc) ;
2793 Row [r].start = (int) (pdest - A [0]) ;
2794 length = Row [r].length ;
2795 for (j = 0 ; j < length ; j++)
2797 c = psrc++ ;
2798 if (COL_IS_ALIVE (Col, c))
2800 pdest = c ;
2801 // pdest++ = c ; // FIXME: pointer
2804 Row [r].length = (int) (pdest - A [Row [r].start]) ;
2805 ASSERT (Row [r].length > 0) ;
2806 if (!NDEBUG)
2808 debug_rows-- ;
2809 } /* NDEBUG */
2812 /* ensure we found all the rows */
2813 ASSERT (debug_rows == 0) ;
2815 /* === Return the new value of pfree ================================ */
2817 return ((int) (pdest - A [0])) ;
2821 /* ========================================================================== */
2822 /* === clear_mark =========================================================== */
2823 /* ========================================================================== */
2826 * Clears the Row [].shared2.mark array, and returns the new tag_mark.
2828 * @param tag_mark new value of tag_mark
2829 * @param max_mark max allowed value of tag_mark
2830 * @param n_row number of rows in A
2831 * @param Row Row [0 ... n_row-1].shared2.mark is set to zero
2832 * @return the new value for tag_mark
2834 private static int clear_mark (int tag_mark, int max_mark, int n_row,
2835 Colamd_Row[] Row)
2837 int r ;
2839 if (tag_mark <= 0 || tag_mark >= max_mark)
2841 for (r = 0 ; r < n_row ; r++)
2843 if (ROW_IS_ALIVE (Row, r))
2845 Row [r].mark = 0 ;
2848 tag_mark = 1 ;
2851 return (tag_mark) ;
2855 /* ========================================================================== */
2856 /* === print_report ========================================================= */
2857 /* ========================================================================== */
2861 * @param method
2862 * @param stats
2864 private static void print_report (String method, int[] stats)
2866 int i1, i2, i3 ;
2868 PRINTF ("\n%s version %d.%d, %s: ", method,
2869 COLAMD_MAIN_VERSION, COLAMD_SUB_VERSION, COLAMD_DATE) ;
2871 if (stats == null)
2873 PRINTF ("No statistics available.\n") ;
2874 return ;
2877 i1 = stats [COLAMD_INFO1] ;
2878 i2 = stats [COLAMD_INFO2] ;
2879 i3 = stats [COLAMD_INFO3] ;
2881 if (stats [COLAMD_STATUS] >= 0)
2883 PRINTF ("OK. ") ;
2885 else
2887 PRINTF ("ERROR. ") ;
2890 switch (stats [COLAMD_STATUS])
2893 case COLAMD_OK_BUT_JUMBLED:
2895 PRINTF("Matrix has unsorted or duplicate row indices.\n") ;
2897 PRINTF("%s: number of duplicate or out-of-order row indices: %d\n",
2898 method, i3) ;
2900 PRINTF("%s: last seen duplicate or out-of-order row index: %d\n",
2901 method, INDEX (i2)) ;
2903 PRINTF("%s: last seen in column: %d",
2904 method, INDEX (i1)) ;
2906 /* no break - fall through to next case instead */
2908 case COLAMD_OK:
2910 PRINTF("\n") ;
2912 PRINTF("%s: number of dense or empty rows ignored: %d\n",
2913 method, stats [COLAMD_DENSE_ROW]) ;
2915 PRINTF("%s: number of dense or empty columns ignored: %d\n",
2916 method, stats [COLAMD_DENSE_COL]) ;
2918 PRINTF("%s: number of garbage collections performed: %d\n",
2919 method, stats [COLAMD_DEFRAG_COUNT]) ;
2920 break ;
2922 case COLAMD_ERROR_A_not_present:
2924 PRINTF("Array A (row indices of matrix) not present.\n") ;
2925 break ;
2927 case COLAMD_ERROR_p_not_present:
2929 PRINTF("Array p (column pointers for matrix) not present.\n") ;
2930 break ;
2932 case COLAMD_ERROR_nrow_negative:
2934 PRINTF("Invalid number of rows (%d).\n", i1) ;
2935 break ;
2937 case COLAMD_ERROR_ncol_negative:
2939 PRINTF("Invalid number of columns (%d).\n", i1) ;
2940 break ;
2942 case COLAMD_ERROR_nnz_negative:
2944 PRINTF("Invalid number of nonzero entries (%d).\n", i1) ;
2945 break ;
2947 case COLAMD_ERROR_p0_nonzero:
2949 PRINTF("Invalid column pointer, p [0] = %d, must be zero.\n", i1);
2950 break ;
2952 case COLAMD_ERROR_A_too_small:
2954 PRINTF("Array A too small.\n") ;
2955 PRINTF(" Need Alen >= %d, but given only Alen = %d.\n",
2956 i1, i2) ;
2957 break ;
2959 case COLAMD_ERROR_col_length_negative:
2961 PRINTF
2962 ("Column %d has a negative number of nonzero entries (%d).\n",
2963 INDEX (i1), i2) ;
2964 break ;
2966 case COLAMD_ERROR_row_index_out_of_bounds:
2968 PRINTF
2969 ("Row index (row %d) out of bounds (%d to %d) in column %d.\n",
2970 INDEX (i2), INDEX (0), INDEX (i3-1), INDEX (i1)) ;
2971 break ;
2973 case COLAMD_ERROR_out_of_memory:
2975 PRINTF("Out of memory.\n") ;
2976 break ;
2978 /* v2.4: internal-error case deleted */
2983 /* ========================================================================== */
2984 /* === colamd debugging routines ============================================ */
2985 /* ========================================================================== */
2988 * At this point, all empty rows and columns are dead. All live columns
2989 * are "clean" (containing no dead rows) and simplicial (no supercolumns
2990 * yet). Rows may contain dead columns, but all live rows contain at
2991 * least one live column.
2993 * @param n_row
2994 * @param n_col
2995 * @param Row
2996 * @param Col
2997 * @param A
2998 * @param n_col2
3000 private static void debug_structures (int n_row, int n_col,
3001 Colamd_Row[] Row, Colamd_Col[] Col, int[] A, int n_col2)
3003 if (!NDEBUG)
3005 /* === Local variables ============================================== */
3007 int i ;
3008 int c ;
3009 int cp ;
3010 int cp_end ;
3011 int len ;
3012 int score ;
3013 int r ;
3014 int rp ;
3015 int rp_end ;
3016 int deg ;
3018 /* === Check A, Row, and Col ======================================== */
3020 for (c = 0 ; c < n_col ; c++)
3022 if (COL_IS_ALIVE (Col, c))
3024 len = Col [c].length ;
3025 score = Col [c].score ;
3026 DEBUG4 ("initial live col %5d %5d %5d\n", c, len, score) ;
3027 ASSERT (len > 0) ;
3028 ASSERT (score >= 0) ;
3029 ASSERT (Col [c].thickness == 1) ;
3030 cp = A [Col [c].start] ;
3031 cp_end = cp + len ;
3032 while (cp < cp_end)
3034 r = cp++ ;
3035 ASSERT (ROW_IS_ALIVE (Row, r)) ;
3038 else
3040 i = Col [c].order ;
3041 ASSERT (i >= n_col2 && i < n_col) ;
3045 for (r = 0 ; r < n_row ; r++)
3047 if (ROW_IS_ALIVE (Row, r))
3049 i = 0 ;
3050 len = Row [r].length ;
3051 deg = Row [r].degree ;
3052 ASSERT (len > 0) ;
3053 ASSERT (deg > 0) ;
3054 rp = A [Row [r].start] ;
3055 rp_end = rp + len ;
3056 while (rp < rp_end)
3058 c = rp++ ;
3059 if (COL_IS_ALIVE (Col, c))
3061 i++ ;
3064 ASSERT (i > 0) ;
3071 /* ========================================================================== */
3072 /* === debug_deg_lists ====================================================== */
3073 /* ========================================================================== */
3076 * Prints the contents of the degree lists. Counts the number of columns
3077 * in the degree list and compares it to the total it should have. Also
3078 * checks the row degrees.
3080 * @param n_row
3081 * @param n_col
3082 * @param Row
3083 * @param Col
3084 * @param head
3085 * @param min_score
3086 * @param should
3087 * @param max_deg
3089 private static void debug_deg_lists (int n_row, int n_col,
3090 Colamd_Row[] Row, Colamd_Col[] Col, int[] head, int min_score,
3091 int should, int max_deg)
3093 if (!NDEBUG)
3095 /* === Local variables ============================================== */
3097 int deg ;
3098 int col ;
3099 int have ;
3100 int row ;
3102 /* === Check the degree lists ======================================= */
3104 if (n_col > 10000 && colamd_debug <= 0)
3106 return ;
3108 have = 0 ;
3109 DEBUG4 ("Degree lists: %d\n", min_score) ;
3110 for (deg = 0 ; deg <= n_col ; deg++)
3112 col = head [deg] ;
3113 if (col == EMPTY)
3115 continue ;
3117 DEBUG4 ("%d:", deg) ;
3118 while (col != EMPTY)
3120 DEBUG4 (" %d", col) ;
3121 have += Col [col].thickness ;
3122 ASSERT (COL_IS_ALIVE (Col, col)) ;
3123 col = Col [col].degree_next ;
3125 DEBUG4 ("\n") ;
3127 DEBUG4 ("should %d have %d\n", should, have) ;
3128 ASSERT (should == have) ;
3130 /* === Check the row degrees ======================================== */
3132 if (n_row > 10000 && colamd_debug <= 0)
3134 return ;
3136 for (row = 0 ; row < n_row ; row++)
3138 if (ROW_IS_ALIVE (Row, row))
3140 ASSERT (Row [row].degree <= max_deg) ;
3147 /* ========================================================================== */
3148 /* === debug_mark =========================================================== */
3149 /* ========================================================================== */
3152 * Ensures that the tag_mark is less that the maximum and also ensures that
3153 * each entry in the mark array is less than the tag mark.
3155 * @param n_row
3156 * @param Row
3157 * @param tag_mark
3158 * @param max_mark
3160 private static void debug_mark (int n_row, Colamd_Row[] Row, int tag_mark,
3161 int max_mark)
3163 if (!NDEBUG)
3165 /* === Local variables ============================================== */
3167 int r ;
3169 /* === Check the Row marks ========================================== */
3171 ASSERT (tag_mark > 0 && tag_mark <= max_mark) ;
3172 if (n_row > 10000 && colamd_debug <= 0)
3174 return ;
3176 for (r = 0 ; r < n_row ; r++)
3178 ASSERT (Row [r].mark < tag_mark) ;
3184 /* ========================================================================== */
3185 /* === debug_matrix ========================================================= */
3186 /* ========================================================================== */
3189 * Prints out the contents of the columns and the rows.
3191 * @param n_row
3192 * @param n_col
3193 * @param Row
3194 * @param Col
3195 * @param A
3197 private static void debug_matrix (int n_row, int n_col,
3198 Colamd_Row[] Row, Colamd_Col[] Col, int[] A)
3200 if (!NDEBUG)
3202 /* === Local variables ============================================== */
3204 int r ;
3205 int c ;
3206 int rp ;
3207 int rp_end ;
3208 int cp ;
3209 int cp_end ;
3211 /* === Dump the rows and columns of the matrix ====================== */
3213 if (colamd_debug < 3)
3215 return ;
3217 DEBUG3 ("DUMP MATRIX:\n") ;
3218 for (r = 0 ; r < n_row ; r++)
3220 DEBUG3 ("Row %d alive? %d\n", r, ROW_IS_ALIVE (Row, r)) ;
3221 if (ROW_IS_DEAD (Row, r))
3223 continue ;
3225 DEBUG3 ("start %d length %d degree %d\n",
3226 Row [r].start, Row [r].length, Row [r].degree) ;
3227 rp = A [Row [r].start] ;
3228 rp_end = rp + Row [r].length ;
3229 while (rp < rp_end)
3231 c = rp++ ;
3232 DEBUG4 (" %d col %d\n", COL_IS_ALIVE (Col, c), c) ;
3236 for (c = 0 ; c < n_col ; c++)
3238 DEBUG3 ("Col %d alive? %d\n", c, COL_IS_ALIVE (Col, c)) ;
3239 if (COL_IS_DEAD (Col, c))
3241 continue ;
3243 DEBUG3 ("start %d length %d shared1 %d shared2 %d\n",
3244 Col [c].start, Col [c].length,
3245 Col [c].thickness, Col [c].score) ;
3246 cp = A [Col [c].start] ;
3247 cp_end = cp + Col [c].length ;
3248 while (cp < cp_end)
3250 r = cp++ ;
3251 DEBUG4 (" %d row %d\n", ROW_IS_ALIVE (Row, r), r) ;
3257 @SuppressWarnings("unused")
3258 private static void colamd_get_debug (String method)
3260 if (!NDEBUG)
3262 File f ;
3263 colamd_debug = 0 ; /* no debug printing */
3264 f = new File("debug") ;
3265 if (!f.exists())
3267 colamd_debug = 0 ;
3269 else
3271 try {
3272 FileReader fr ;
3273 fr = new FileReader(f) ;
3274 BufferedReader br ;
3275 br = new BufferedReader(fr) ;
3276 colamd_debug = Integer.valueOf( br.readLine() ) ;
3277 br.close() ;
3278 fr.close() ;
3279 } catch (IOException e) {
3280 PRINTF ("%s: AMD_debug_init, " +
3281 "error reading debug.amd file", method) ;
3284 DEBUG0 ("%s: debug version, D = %d (THIS WILL BE SLOW!)\n",
3285 method, colamd_debug) ;