Eigen_Colamd.h
Go to the documentation of this file.
1 // // This file is part of Eigen, a lightweight C++ template library
2 // for linear algebra.
3 //
4 // Copyright (C) 2012 Desire Nuentsa Wakam <desire.nuentsa_wakam@inria.fr>
5 //
6 // This Source Code Form is subject to the terms of the Mozilla
7 // Public License v. 2.0. If a copy of the MPL was not distributed
8 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
9 
10 // This file is modified from the colamd/symamd library. The copyright is below
11 
12 // The authors of the code itself are Stefan I. Larimore and Timothy A.
13 // Davis (davis@cise.ufl.edu), University of Florida. The algorithm was
14 // developed in collaboration with John Gilbert, Xerox PARC, and Esmond
15 // Ng, Oak Ridge National Laboratory.
16 //
17 // Date:
18 //
19 // September 8, 2003. Version 2.3.
20 //
21 // Acknowledgements:
22 //
23 // This work was supported by the National Science Foundation, under
24 // grants DMS-9504974 and DMS-9803599.
25 //
26 // Notice:
27 //
28 // Copyright (c) 1998-2003 by the University of Florida.
29 // All Rights Reserved.
30 //
31 // THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY
32 // EXPRESSED OR IMPLIED. ANY USE IS AT YOUR OWN RISK.
33 //
34 // Permission is hereby granted to use, copy, modify, and/or distribute
35 // this program, provided that the Copyright, this License, and the
36 // Availability of the original version is retained on all copies and made
37 // accessible to the end-user of any code or package that includes COLAMD
38 // or any modified version of COLAMD.
39 //
40 // Availability:
41 //
42 // The colamd/symamd library is available at
43 //
44 // http://www.suitesparse.com
45 
46 #ifndef EIGEN_COLAMD_H
47 #define EIGEN_COLAMD_H
48 
49 namespace internal {
50 
51 namespace Colamd {
52 
53 /* Ensure that debugging is turned off: */
54 #ifndef COLAMD_NDEBUG
55 #define COLAMD_NDEBUG
56 #endif /* NDEBUG */
57 
58 /* ========================================================================== */
59 /* === Knob and statistics definitions ====================================== */
60 /* ========================================================================== */
61 
62 /* size of the knobs [ ] array. Only knobs [0..1] are currently used. */
63 const int NKnobs = 20;
64 
65 /* number of output statistics. Only stats [0..6] are currently used. */
66 const int NStats = 20;
67 
68 /* Indices into knobs and stats array. */
70  /* knobs [0] and stats [0]: dense row knob and output statistic. */
71  DenseRow = 0,
72 
73  /* knobs [1] and stats [1]: dense column knob and output statistic. */
74  DenseCol = 1,
75 
76  /* stats [2]: memory defragmentation count output statistic */
78 
79  /* stats [3]: colamd status: zero OK, > 0 warning or notice, < 0 error */
80  Status = 3,
81 
82  /* stats [4..6]: error info, or info on jumbled columns */
83  Info1 = 4,
84  Info2 = 5,
85  Info3 = 6
86 };
87 
88 /* error codes returned in stats [3]: */
89 enum Status {
90  Ok = 0,
102  ErrorInternalError = -999
103 };
104 /* ========================================================================== */
105 /* === Definitions ========================================================== */
106 /* ========================================================================== */
107 
108 template <typename IndexType>
109 IndexType ones_complement(const IndexType r) {
110  return (-(r)-1);
111 }
112 
113 /* -------------------------------------------------------------------------- */
114 const int Empty = -1;
115 
116 /* Row and column status */
117 enum RowColumnStatus { Alive = 0, Dead = -1 };
118 
119 /* Column status */
121 
122 /* ========================================================================== */
123 /* === Colamd reporting mechanism =========================================== */
124 /* ========================================================================== */
125 
126 // == Row and Column structures ==
127 template <typename IndexType>
128 struct ColStructure {
129  IndexType start; /* index for A of first row in this column, or Dead */
130  /* if column is dead */
131  IndexType length; /* number of rows in this column */
132  union {
133  IndexType thickness; /* number of original columns represented by this */
134  /* col, if the column is alive */
135  IndexType parent; /* parent in parent tree super-column structure, if */
136  /* the column is dead */
138  union {
139  IndexType score; /* the score used to maintain heap, if col is alive */
140  IndexType order; /* pivot ordering of this column, if col is dead */
142  union {
143  IndexType headhash; /* head of a hash bucket, if col is at the head of */
144  /* a degree list */
145  IndexType hash; /* hash value, if col is not in a degree list */
146  IndexType prev; /* previous column in degree list, if col is in a */
147  /* degree list (but not at the head of a degree list) */
149  union {
150  IndexType degree_next; /* next column, if col is in a degree list */
151  IndexType hash_next; /* next column, if col is in a hash list */
153 
154  inline bool is_dead() const { return start < Alive; }
155 
156  inline bool is_alive() const { return start >= Alive; }
157 
158  inline bool is_dead_principal() const { return start == DeadPrincipal; }
159 
160  inline void kill_principal() { start = DeadPrincipal; }
161 
163 };
164 
165 template <typename IndexType>
166 struct RowStructure {
167  IndexType start; /* index for A of first col in this row */
168  IndexType length; /* number of principal columns in this row */
169  union {
170  IndexType degree; /* number of principal & non-principal columns in row */
171  IndexType p; /* used as a row pointer in init_rows_cols () */
173  union {
174  IndexType mark; /* for computing set differences and marking dead rows*/
175  IndexType first_column; /* first column in row (used in garbage collection) */
177 
178  inline bool is_dead() const { return shared2.mark < Alive; }
179 
180  inline bool is_alive() const { return shared2.mark >= Alive; }
181 
182  inline void kill() { shared2.mark = Dead; }
183 };
184 
185 /* ========================================================================== */
186 /* === Colamd recommended memory size ======================================= */
187 /* ========================================================================== */
188 
189 /*
190  The recommended length Alen of the array A passed to colamd is given by
191  the COLAMD_RECOMMENDED (nnz, n_row, n_col) macro. It returns -1 if any
192  argument is negative. 2*nnz space is required for the row and column
193  indices of the matrix. colamd_c (n_col) + colamd_r (n_row) space is
194  required for the Col and Row arrays, respectively, which are internal to
195  colamd. An additional n_col space is the minimal amount of "elbow room",
196  and nnz/5 more space is recommended for run time efficiency.
197 
198  This macro is not needed when using symamd.
199 
200  Explicit typecast to IndexType added Sept. 23, 2002, COLAMD version 2.2, to avoid
201  gcc -pedantic warning messages.
202 */
203 template <typename IndexType>
204 inline IndexType colamd_c(IndexType n_col) {
205  return IndexType(((n_col) + 1) * sizeof(ColStructure<IndexType>) / sizeof(IndexType));
206 }
207 
208 template <typename IndexType>
209 inline IndexType colamd_r(IndexType n_row) {
210  return IndexType(((n_row) + 1) * sizeof(RowStructure<IndexType>) / sizeof(IndexType));
211 }
212 
213 // Prototypes of non-user callable routines
214 template <typename IndexType>
215 static IndexType init_rows_cols(IndexType n_row, IndexType n_col, RowStructure<IndexType> Row[],
216  ColStructure<IndexType> col[], IndexType A[], IndexType p[], IndexType stats[NStats]);
217 
218 template <typename IndexType>
219 static void init_scoring(IndexType n_row, IndexType n_col, RowStructure<IndexType> Row[], ColStructure<IndexType> Col[],
220  IndexType A[], IndexType head[], double knobs[NKnobs], IndexType *p_n_row2,
221  IndexType *p_n_col2, IndexType *p_max_deg);
222 
223 template <typename IndexType>
224 static IndexType find_ordering(IndexType n_row, IndexType n_col, IndexType Alen, RowStructure<IndexType> Row[],
225  ColStructure<IndexType> Col[], IndexType A[], IndexType head[], IndexType n_col2,
226  IndexType max_deg, IndexType pfree);
227 
228 template <typename IndexType>
229 static void order_children(IndexType n_col, ColStructure<IndexType> Col[], IndexType p[]);
230 
231 template <typename IndexType>
232 static void detect_super_cols(ColStructure<IndexType> Col[], IndexType A[], IndexType head[], IndexType row_start,
233  IndexType row_length);
234 
235 template <typename IndexType>
236 static IndexType garbage_collection(IndexType n_row, IndexType n_col, RowStructure<IndexType> Row[],
237  ColStructure<IndexType> Col[], IndexType A[], IndexType *pfree);
238 
239 template <typename IndexType>
240 static inline IndexType clear_mark(IndexType n_row, RowStructure<IndexType> Row[]);
241 
242 /* === No debugging ========================================================= */
243 
244 #define COLAMD_DEBUG0(params) ;
245 #define COLAMD_DEBUG1(params) ;
246 #define COLAMD_DEBUG2(params) ;
247 #define COLAMD_DEBUG3(params) ;
248 #define COLAMD_DEBUG4(params) ;
249 
250 #define COLAMD_ASSERT(expression) ((void)0)
251 
266 template <typename IndexType>
267 inline IndexType recommended(IndexType nnz, IndexType n_row, IndexType n_col) {
268  if ((nnz) < 0 || (n_row) < 0 || (n_col) < 0)
269  return (-1);
270  else
271  return (2 * (nnz) + colamd_c(n_col) + colamd_r(n_row) + (n_col) + ((nnz) / 5));
272 }
273 
295 static inline void set_defaults(double knobs[NKnobs]) {
296  /* === Local variables ================================================== */
297 
298  int i;
299 
300  if (!knobs) {
301  return; /* no knobs to initialize */
302  }
303  for (i = 0; i < NKnobs; i++) {
304  knobs[i] = 0;
305  }
306  knobs[Colamd::DenseRow] = 0.5; /* ignore rows over 50% dense */
307  knobs[Colamd::DenseCol] = 0.5; /* ignore columns over 50% dense */
308 }
309 
327 template <typename IndexType>
328 static bool compute_ordering(IndexType n_row, IndexType n_col, IndexType Alen, IndexType *A, IndexType *p,
329  double knobs[NKnobs], IndexType stats[NStats]) {
330  /* === Local variables ================================================== */
331 
332  IndexType i; /* loop index */
333  IndexType nnz; /* nonzeros in A */
334  IndexType Row_size; /* size of Row [], in integers */
335  IndexType Col_size; /* size of Col [], in integers */
336  IndexType need; /* minimum required length of A */
337  Colamd::RowStructure<IndexType> *Row; /* pointer into A of Row [0..n_row] array */
338  Colamd::ColStructure<IndexType> *Col; /* pointer into A of Col [0..n_col] array */
339  IndexType n_col2; /* number of non-dense, non-empty columns */
340  IndexType n_row2; /* number of non-dense, non-empty rows */
341  IndexType ngarbage; /* number of garbage collections performed */
342  IndexType max_deg; /* maximum row degree */
343  double default_knobs[NKnobs]; /* default knobs array */
344 
345  /* === Check the input arguments ======================================== */
346 
347  if (!stats) {
348  COLAMD_DEBUG0(("colamd: stats not present\n"));
349  return (false);
350  }
351  for (i = 0; i < NStats; i++) {
352  stats[i] = 0;
353  }
354  stats[Colamd::Status] = Colamd::Ok;
355  stats[Colamd::Info1] = -1;
356  stats[Colamd::Info2] = -1;
357 
358  if (!A) /* A is not present */
359  {
361  COLAMD_DEBUG0(("colamd: A not present\n"));
362  return (false);
363  }
364 
365  if (!p) /* p is not present */
366  {
368  COLAMD_DEBUG0(("colamd: p not present\n"));
369  return (false);
370  }
371 
372  if (n_row < 0) /* n_row must be >= 0 */
373  {
375  stats[Colamd::Info1] = n_row;
376  COLAMD_DEBUG0(("colamd: nrow negative %d\n", n_row));
377  return (false);
378  }
379 
380  if (n_col < 0) /* n_col must be >= 0 */
381  {
383  stats[Colamd::Info1] = n_col;
384  COLAMD_DEBUG0(("colamd: ncol negative %d\n", n_col));
385  return (false);
386  }
387 
388  nnz = p[n_col];
389  if (nnz < 0) /* nnz must be >= 0 */
390  {
392  stats[Colamd::Info1] = nnz;
393  COLAMD_DEBUG0(("colamd: number of entries negative %d\n", nnz));
394  return (false);
395  }
396 
397  if (p[0] != 0) {
399  stats[Colamd::Info1] = p[0];
400  COLAMD_DEBUG0(("colamd: p[0] not zero %d\n", p[0]));
401  return (false);
402  }
403 
404  /* === If no knobs, set default knobs =================================== */
405 
406  if (!knobs) {
407  set_defaults(default_knobs);
408  knobs = default_knobs;
409  }
410 
411  /* === Allocate the Row and Col arrays from array A ===================== */
412 
413  Col_size = colamd_c(n_col);
414  Row_size = colamd_r(n_row);
415  need = 2 * nnz + n_col + Col_size + Row_size;
416 
417  if (need > Alen) {
418  /* not enough space in array A to perform the ordering */
420  stats[Colamd::Info1] = need;
421  stats[Colamd::Info2] = Alen;
422  COLAMD_DEBUG0(("colamd: Need Alen >= %d, given only Alen = %d\n", need, Alen));
423  return (false);
424  }
425 
426  Alen -= Col_size + Row_size;
427  Col = (ColStructure<IndexType> *)&A[Alen];
428  Row = (RowStructure<IndexType> *)&A[Alen + Col_size];
429 
430  /* === Construct the row and column data structures ===================== */
431 
432  if (!Colamd::init_rows_cols(n_row, n_col, Row, Col, A, p, stats)) {
433  /* input matrix is invalid */
434  COLAMD_DEBUG0(("colamd: Matrix invalid\n"));
435  return (false);
436  }
437 
438  /* === Initialize scores, kill dense rows/columns ======================= */
439 
440  Colamd::init_scoring(n_row, n_col, Row, Col, A, p, knobs, &n_row2, &n_col2, &max_deg);
441 
442  /* === Order the supercolumns =========================================== */
443 
444  ngarbage = Colamd::find_ordering(n_row, n_col, Alen, Row, Col, A, p, n_col2, max_deg, 2 * nnz);
445 
446  /* === Order the non-principal columns ================================== */
447 
448  Colamd::order_children(n_col, Col, p);
449 
450  /* === Return statistics in stats ======================================= */
451 
452  stats[Colamd::DenseRow] = n_row - n_row2;
453  stats[Colamd::DenseCol] = n_col - n_col2;
454  stats[Colamd::DefragCount] = ngarbage;
455  COLAMD_DEBUG0(("colamd: done.\n"));
456  return (true);
457 }
458 
459 /* ========================================================================== */
460 /* === NON-USER-CALLABLE ROUTINES: ========================================== */
461 /* ========================================================================== */
462 
463 /* There are no user-callable routines beyond this point in the file */
464 
465 /* ========================================================================== */
466 /* === init_rows_cols ======================================================= */
467 /* ========================================================================== */
468 
469 /*
470  Takes the column form of the matrix in A and creates the row form of the
471  matrix. Also, row and column attributes are stored in the Col and Row
472  structs. If the columns are un-sorted or contain duplicate row indices,
473  this routine will also sort and remove duplicate row indices from the
474  column form of the matrix. Returns false if the matrix is invalid,
475  true otherwise. Not user-callable.
476 */
477 template <typename IndexType>
478 static IndexType init_rows_cols /* returns true if OK, or false otherwise */
479  (
480  /* === Parameters ======================================================= */
481 
482  IndexType n_row, /* number of rows of A */
483  IndexType n_col, /* number of columns of A */
484  RowStructure<IndexType> Row[], /* of size n_row+1 */
485  ColStructure<IndexType> Col[], /* of size n_col+1 */
486  IndexType A[], /* row indices of A, of size Alen */
487  IndexType p[], /* pointers to columns in A, of size n_col+1 */
488  IndexType stats[NStats] /* colamd statistics */
489  ) {
490  /* === Local variables ================================================== */
491 
492  IndexType col; /* a column index */
493  IndexType row; /* a row index */
494  IndexType *cp; /* a column pointer */
495  IndexType *cp_end; /* a pointer to the end of a column */
496  IndexType *rp; /* a row pointer */
497  IndexType *rp_end; /* a pointer to the end of a row */
498  IndexType last_row; /* previous row */
499 
500  /* === Initialize columns, and check column pointers ==================== */
501 
502  for (col = 0; col < n_col; col++) {
503  Col[col].start = p[col];
504  Col[col].length = p[col + 1] - p[col];
505 
506  if ((Col[col].length) < 0) // extra parentheses to work-around gcc bug 10200
507  {
508  /* column pointers must be non-decreasing */
510  stats[Colamd::Info1] = col;
511  stats[Colamd::Info2] = Col[col].length;
512  COLAMD_DEBUG0(("colamd: col %d length %d < 0\n", col, Col[col].length));
513  return (false);
514  }
515 
516  Col[col].shared1.thickness = 1;
517  Col[col].shared2.score = 0;
518  Col[col].shared3.prev = Empty;
519  Col[col].shared4.degree_next = Empty;
520  }
521 
522  /* p [0..n_col] no longer needed, used as "head" in subsequent routines */
523 
524  /* === Scan columns, compute row degrees, and check row indices ========= */
525 
526  stats[Info3] = 0; /* number of duplicate or unsorted row indices*/
527 
528  for (row = 0; row < n_row; row++) {
529  Row[row].length = 0;
530  Row[row].shared2.mark = -1;
531  }
532 
533  for (col = 0; col < n_col; col++) {
534  last_row = -1;
535 
536  cp = &A[p[col]];
537  cp_end = &A[p[col + 1]];
538 
539  while (cp < cp_end) {
540  row = *cp++;
541 
542  /* make sure row indices within range */
543  if (row < 0 || row >= n_row) {
545  stats[Colamd::Info1] = col;
546  stats[Colamd::Info2] = row;
547  stats[Colamd::Info3] = n_row;
548  COLAMD_DEBUG0(("colamd: row %d col %d out of bounds\n", row, col));
549  return (false);
550  }
551 
552  if (row <= last_row || Row[row].shared2.mark == col) {
553  /* row index are unsorted or repeated (or both), thus col */
554  /* is jumbled. This is a notice, not an error condition. */
556  stats[Colamd::Info1] = col;
557  stats[Colamd::Info2] = row;
558  (stats[Colamd::Info3])++;
559  COLAMD_DEBUG1(("colamd: row %d col %d unsorted/duplicate\n", row, col));
560  }
561 
562  if (Row[row].shared2.mark != col) {
563  Row[row].length++;
564  } else {
565  /* this is a repeated entry in the column, */
566  /* it will be removed */
567  Col[col].length--;
568  }
569 
570  /* mark the row as having been seen in this column */
571  Row[row].shared2.mark = col;
572 
573  last_row = row;
574  }
575  }
576 
577  /* === Compute row pointers ============================================= */
578 
579  /* row form of the matrix starts directly after the column */
580  /* form of matrix in A */
581  Row[0].start = p[n_col];
582  Row[0].shared1.p = Row[0].start;
583  Row[0].shared2.mark = -1;
584  for (row = 1; row < n_row; row++) {
585  Row[row].start = Row[row - 1].start + Row[row - 1].length;
586  Row[row].shared1.p = Row[row].start;
587  Row[row].shared2.mark = -1;
588  }
589 
590  /* === Create row form ================================================== */
591 
592  if (stats[Status] == OkButJumbled) {
593  /* if cols jumbled, watch for repeated row indices */
594  for (col = 0; col < n_col; col++) {
595  cp = &A[p[col]];
596  cp_end = &A[p[col + 1]];
597  while (cp < cp_end) {
598  row = *cp++;
599  if (Row[row].shared2.mark != col) {
600  A[(Row[row].shared1.p)++] = col;
601  Row[row].shared2.mark = col;
602  }
603  }
604  }
605  } else {
606  /* if cols not jumbled, we don't need the mark (this is faster) */
607  for (col = 0; col < n_col; col++) {
608  cp = &A[p[col]];
609  cp_end = &A[p[col + 1]];
610  while (cp < cp_end) {
611  A[(Row[*cp++].shared1.p)++] = col;
612  }
613  }
614  }
615 
616  /* === Clear the row marks and set row degrees ========================== */
617 
618  for (row = 0; row < n_row; row++) {
619  Row[row].shared2.mark = 0;
620  Row[row].shared1.degree = Row[row].length;
621  }
622 
623  /* === See if we need to re-create columns ============================== */
624 
625  if (stats[Status] == OkButJumbled) {
626  COLAMD_DEBUG0(("colamd: reconstructing column form, matrix jumbled\n"));
627 
628  /* === Compute col pointers ========================================= */
629 
630  /* col form of the matrix starts at A [0]. */
631  /* Note, we may have a gap between the col form and the row */
632  /* form if there were duplicate entries, if so, it will be */
633  /* removed upon the first garbage collection */
634  Col[0].start = 0;
635  p[0] = Col[0].start;
636  for (col = 1; col < n_col; col++) {
637  /* note that the lengths here are for pruned columns, i.e. */
638  /* no duplicate row indices will exist for these columns */
639  Col[col].start = Col[col - 1].start + Col[col - 1].length;
640  p[col] = Col[col].start;
641  }
642 
643  /* === Re-create col form =========================================== */
644 
645  for (row = 0; row < n_row; row++) {
646  rp = &A[Row[row].start];
647  rp_end = rp + Row[row].length;
648  while (rp < rp_end) {
649  A[(p[*rp++])++] = row;
650  }
651  }
652  }
653 
654  /* === Done. Matrix is not (or no longer) jumbled ====================== */
655 
656  return (true);
657 }
658 
659 /* ========================================================================== */
660 /* === init_scoring ========================================================= */
661 /* ========================================================================== */
662 
663 /*
664  Kills dense or empty columns and rows, calculates an initial score for
665  each column, and places all columns in the degree lists. Not user-callable.
666 */
667 template <typename IndexType>
668 static void init_scoring(
669  /* === Parameters ======================================================= */
670 
671  IndexType n_row, /* number of rows of A */
672  IndexType n_col, /* number of columns of A */
673  RowStructure<IndexType> Row[], /* of size n_row+1 */
674  ColStructure<IndexType> Col[], /* of size n_col+1 */
675  IndexType A[], /* column form and row form of A */
676  IndexType head[], /* of size n_col+1 */
677  double knobs[NKnobs], /* parameters */
678  IndexType *p_n_row2, /* number of non-dense, non-empty rows */
679  IndexType *p_n_col2, /* number of non-dense, non-empty columns */
680  IndexType *p_max_deg /* maximum row degree */
681 ) {
682  /* === Local variables ================================================== */
683 
684  IndexType c; /* a column index */
685  IndexType r, row; /* a row index */
686  IndexType *cp; /* a column pointer */
687  IndexType deg; /* degree of a row or column */
688  IndexType *cp_end; /* a pointer to the end of a column */
689  IndexType *new_cp; /* new column pointer */
690  IndexType col_length; /* length of pruned column */
691  IndexType score; /* current column score */
692  IndexType n_col2; /* number of non-dense, non-empty columns */
693  IndexType n_row2; /* number of non-dense, non-empty rows */
694  IndexType dense_row_count; /* remove rows with more entries than this */
695  IndexType dense_col_count; /* remove cols with more entries than this */
696  IndexType min_score; /* smallest column score */
697  IndexType max_deg; /* maximum row degree */
698  IndexType next_col; /* Used to add to degree list.*/
699 
700  /* === Extract knobs ==================================================== */
701 
702  dense_row_count = numext::maxi(IndexType(0), numext::mini(IndexType(knobs[Colamd::DenseRow] * n_col), n_col));
703  dense_col_count = numext::maxi(IndexType(0), numext::mini(IndexType(knobs[Colamd::DenseCol] * n_row), n_row));
704  COLAMD_DEBUG1(("colamd: densecount: %d %d\n", dense_row_count, dense_col_count));
705  max_deg = 0;
706  n_col2 = n_col;
707  n_row2 = n_row;
708 
709  /* === Kill empty columns =============================================== */
710 
711  /* Put the empty columns at the end in their natural order, so that LU */
712  /* factorization can proceed as far as possible. */
713  for (c = n_col - 1; c >= 0; c--) {
714  deg = Col[c].length;
715  if (deg == 0) {
716  /* this is a empty column, kill and order it last */
717  Col[c].shared2.order = --n_col2;
718  Col[c].kill_principal();
719  }
720  }
721  COLAMD_DEBUG1(("colamd: null columns killed: %d\n", n_col - n_col2));
722 
723  /* === Kill dense columns =============================================== */
724 
725  /* Put the dense columns at the end, in their natural order */
726  for (c = n_col - 1; c >= 0; c--) {
727  /* skip any dead columns */
728  if (Col[c].is_dead()) {
729  continue;
730  }
731  deg = Col[c].length;
732  if (deg > dense_col_count) {
733  /* this is a dense column, kill and order it last */
734  Col[c].shared2.order = --n_col2;
735  /* decrement the row degrees */
736  cp = &A[Col[c].start];
737  cp_end = cp + Col[c].length;
738  while (cp < cp_end) {
739  Row[*cp++].shared1.degree--;
740  }
741  Col[c].kill_principal();
742  }
743  }
744  COLAMD_DEBUG1(("colamd: Dense and null columns killed: %d\n", n_col - n_col2));
745 
746  /* === Kill dense and empty rows ======================================== */
747 
748  for (r = 0; r < n_row; r++) {
749  deg = Row[r].shared1.degree;
750  COLAMD_ASSERT(deg >= 0 && deg <= n_col);
751  if (deg > dense_row_count || deg == 0) {
752  /* kill a dense or empty row */
753  Row[r].kill();
754  --n_row2;
755  } else {
756  /* keep track of max degree of remaining rows */
757  max_deg = numext::maxi(max_deg, deg);
758  }
759  }
760  COLAMD_DEBUG1(("colamd: Dense and null rows killed: %d\n", n_row - n_row2));
761 
762  /* === Compute initial column scores ==================================== */
763 
764  /* At this point the row degrees are accurate. They reflect the number */
765  /* of "live" (non-dense) columns in each row. No empty rows exist. */
766  /* Some "live" columns may contain only dead rows, however. These are */
767  /* pruned in the code below. */
768 
769  /* now find the initial matlab score for each column */
770  for (c = n_col - 1; c >= 0; c--) {
771  /* skip dead column */
772  if (Col[c].is_dead()) {
773  continue;
774  }
775  score = 0;
776  cp = &A[Col[c].start];
777  new_cp = cp;
778  cp_end = cp + Col[c].length;
779  while (cp < cp_end) {
780  /* get a row */
781  row = *cp++;
782  /* skip if dead */
783  if (Row[row].is_dead()) {
784  continue;
785  }
786  /* compact the column */
787  *new_cp++ = row;
788  /* add row's external degree */
789  score += Row[row].shared1.degree - 1;
790  /* guard against integer overflow */
791  score = numext::mini(score, n_col);
792  }
793  /* determine pruned column length */
794  col_length = (IndexType)(new_cp - &A[Col[c].start]);
795  if (col_length == 0) {
796  /* a newly-made null column (all rows in this col are "dense" */
797  /* and have already been killed) */
798  COLAMD_DEBUG2(("Newly null killed: %d\n", c));
799  Col[c].shared2.order = --n_col2;
800  Col[c].kill_principal();
801  } else {
802  /* set column length and set score */
803  COLAMD_ASSERT(score >= 0);
804  COLAMD_ASSERT(score <= n_col);
805  Col[c].length = col_length;
806  Col[c].shared2.score = score;
807  }
808  }
809  COLAMD_DEBUG1(("colamd: Dense, null, and newly-null columns killed: %d\n", n_col - n_col2));
810 
811  /* At this point, all empty rows and columns are dead. All live columns */
812  /* are "clean" (containing no dead rows) and simplicial (no supercolumns */
813  /* yet). Rows may contain dead columns, but all live rows contain at */
814  /* least one live column. */
815 
816  /* === Initialize degree lists ========================================== */
817 
818  /* clear the hash buckets */
819  for (c = 0; c <= n_col; c++) {
820  head[c] = Empty;
821  }
822  min_score = n_col;
823  /* place in reverse order, so low column indices are at the front */
824  /* of the lists. This is to encourage natural tie-breaking */
825  for (c = n_col - 1; c >= 0; c--) {
826  /* only add principal columns to degree lists */
827  if (Col[c].is_alive()) {
828  COLAMD_DEBUG4(("place %d score %d minscore %d ncol %d\n", c, Col[c].shared2.score, min_score, n_col));
829 
830  /* === Add columns score to DList =============================== */
831 
832  score = Col[c].shared2.score;
833 
834  COLAMD_ASSERT(min_score >= 0);
835  COLAMD_ASSERT(min_score <= n_col);
836  COLAMD_ASSERT(score >= 0);
837  COLAMD_ASSERT(score <= n_col);
838  COLAMD_ASSERT(head[score] >= Empty);
839 
840  /* now add this column to dList at proper score location */
841  next_col = head[score];
842  Col[c].shared3.prev = Empty;
843  Col[c].shared4.degree_next = next_col;
844 
845  /* if there already was a column with the same score, set its */
846  /* previous pointer to this new column */
847  if (next_col != Empty) {
848  Col[next_col].shared3.prev = c;
849  }
850  head[score] = c;
851 
852  /* see if this score is less than current min */
853  min_score = numext::mini(min_score, score);
854  }
855  }
856 
857  /* === Return number of remaining columns, and max row degree =========== */
858 
859  *p_n_col2 = n_col2;
860  *p_n_row2 = n_row2;
861  *p_max_deg = max_deg;
862 }
863 
864 /* ========================================================================== */
865 /* === find_ordering ======================================================== */
866 /* ========================================================================== */
867 
868 /*
869  Order the principal columns of the supercolumn form of the matrix
870  (no supercolumns on input). Uses a minimum approximate column minimum
871  degree ordering method. Not user-callable.
872 */
873 template <typename IndexType>
874 static IndexType find_ordering /* return the number of garbage collections */
875  (
876  /* === Parameters ======================================================= */
877 
878  IndexType n_row, /* number of rows of A */
879  IndexType n_col, /* number of columns of A */
880  IndexType Alen, /* size of A, 2*nnz + n_col or larger */
881  RowStructure<IndexType> Row[], /* of size n_row+1 */
882  ColStructure<IndexType> Col[], /* of size n_col+1 */
883  IndexType A[], /* column form and row form of A */
884  IndexType head[], /* of size n_col+1 */
885  IndexType n_col2, /* Remaining columns to order */
886  IndexType max_deg, /* Maximum row degree */
887  IndexType pfree /* index of first free slot (2*nnz on entry) */
888  ) {
889  /* === Local variables ================================================== */
890 
891  IndexType k; /* current pivot ordering step */
892  IndexType pivot_col; /* current pivot column */
893  IndexType *cp; /* a column pointer */
894  IndexType *rp; /* a row pointer */
895  IndexType pivot_row; /* current pivot row */
896  IndexType *new_cp; /* modified column pointer */
897  IndexType *new_rp; /* modified row pointer */
898  IndexType pivot_row_start; /* pointer to start of pivot row */
899  IndexType pivot_row_degree; /* number of columns in pivot row */
900  IndexType pivot_row_length; /* number of supercolumns in pivot row */
901  IndexType pivot_col_score; /* score of pivot column */
902  IndexType needed_memory; /* free space needed for pivot row */
903  IndexType *cp_end; /* pointer to the end of a column */
904  IndexType *rp_end; /* pointer to the end of a row */
905  IndexType row; /* a row index */
906  IndexType col; /* a column index */
907  IndexType max_score; /* maximum possible score */
908  IndexType cur_score; /* score of current column */
909  unsigned int hash; /* hash value for supernode detection */
910  IndexType head_column; /* head of hash bucket */
911  IndexType first_col; /* first column in hash bucket */
912  IndexType tag_mark; /* marker value for mark array */
913  IndexType row_mark; /* Row [row].shared2.mark */
914  IndexType set_difference; /* set difference size of row with pivot row */
915  IndexType min_score; /* smallest column score */
916  IndexType col_thickness; /* "thickness" (no. of columns in a supercol) */
917  IndexType max_mark; /* maximum value of tag_mark */
918  IndexType pivot_col_thickness; /* number of columns represented by pivot col */
919  IndexType prev_col; /* Used by Dlist operations. */
920  IndexType next_col; /* Used by Dlist operations. */
921  IndexType ngarbage; /* number of garbage collections performed */
922 
923  /* === Initialization and clear mark ==================================== */
924 
925  max_mark = INT_MAX - n_col; /* INT_MAX defined in <limits.h> */
926  tag_mark = Colamd::clear_mark(n_row, Row);
927  min_score = 0;
928  ngarbage = 0;
929  COLAMD_DEBUG1(("colamd: Ordering, n_col2=%d\n", n_col2));
930 
931  /* === Order the columns ================================================ */
932 
933  for (k = 0; k < n_col2; /* 'k' is incremented below */) {
934  /* === Select pivot column, and order it ============================ */
935 
936  /* make sure degree list isn't empty */
937  COLAMD_ASSERT(min_score >= 0);
938  COLAMD_ASSERT(min_score <= n_col);
939  COLAMD_ASSERT(head[min_score] >= Empty);
940 
941  /* get pivot column from head of minimum degree list */
942  while (min_score < n_col && head[min_score] == Empty) {
943  min_score++;
944  }
945  pivot_col = head[min_score];
946  COLAMD_ASSERT(pivot_col >= 0 && pivot_col <= n_col);
947  next_col = Col[pivot_col].shared4.degree_next;
948  head[min_score] = next_col;
949  if (next_col != Empty) {
950  Col[next_col].shared3.prev = Empty;
951  }
952 
953  COLAMD_ASSERT(Col[pivot_col].is_alive());
954  COLAMD_DEBUG3(("Pivot col: %d\n", pivot_col));
955 
956  /* remember score for defrag check */
957  pivot_col_score = Col[pivot_col].shared2.score;
958 
959  /* the pivot column is the kth column in the pivot order */
960  Col[pivot_col].shared2.order = k;
961 
962  /* increment order count by column thickness */
963  pivot_col_thickness = Col[pivot_col].shared1.thickness;
964  k += pivot_col_thickness;
965  COLAMD_ASSERT(pivot_col_thickness > 0);
966 
967  /* === Garbage_collection, if necessary ============================= */
968 
969  needed_memory = numext::mini(pivot_col_score, n_col - k);
970  if (pfree + needed_memory >= Alen) {
971  pfree = Colamd::garbage_collection(n_row, n_col, Row, Col, A, &A[pfree]);
972  ngarbage++;
973  /* after garbage collection we will have enough */
974  COLAMD_ASSERT(pfree + needed_memory < Alen);
975  /* garbage collection has wiped out the Row[].shared2.mark array */
976  tag_mark = Colamd::clear_mark(n_row, Row);
977  }
978 
979  /* === Compute pivot row pattern ==================================== */
980 
981  /* get starting location for this new merged row */
982  pivot_row_start = pfree;
983 
984  /* initialize new row counts to zero */
985  pivot_row_degree = 0;
986 
987  /* tag pivot column as having been visited so it isn't included */
988  /* in merged pivot row */
989  Col[pivot_col].shared1.thickness = -pivot_col_thickness;
990 
991  /* pivot row is the union of all rows in the pivot column pattern */
992  cp = &A[Col[pivot_col].start];
993  cp_end = cp + Col[pivot_col].length;
994  while (cp < cp_end) {
995  /* get a row */
996  row = *cp++;
997  COLAMD_DEBUG4(("Pivot col pattern %d %d\n", Row[row].is_alive(), row));
998  /* skip if row is dead */
999  if (Row[row].is_dead()) {
1000  continue;
1001  }
1002  rp = &A[Row[row].start];
1003  rp_end = rp + Row[row].length;
1004  while (rp < rp_end) {
1005  /* get a column */
1006  col = *rp++;
1007  /* add the column, if alive and untagged */
1008  col_thickness = Col[col].shared1.thickness;
1009  if (col_thickness > 0 && Col[col].is_alive()) {
1010  /* tag column in pivot row */
1011  Col[col].shared1.thickness = -col_thickness;
1012  COLAMD_ASSERT(pfree < Alen);
1013  /* place column in pivot row */
1014  A[pfree++] = col;
1015  pivot_row_degree += col_thickness;
1016  }
1017  }
1018  }
1019 
1020  /* clear tag on pivot column */
1021  Col[pivot_col].shared1.thickness = pivot_col_thickness;
1022  max_deg = numext::maxi(max_deg, pivot_row_degree);
1023 
1024  /* === Kill all rows used to construct pivot row ==================== */
1025 
1026  /* also kill pivot row, temporarily */
1027  cp = &A[Col[pivot_col].start];
1028  cp_end = cp + Col[pivot_col].length;
1029  while (cp < cp_end) {
1030  /* may be killing an already dead row */
1031  row = *cp++;
1032  COLAMD_DEBUG3(("Kill row in pivot col: %d\n", row));
1033  Row[row].kill();
1034  }
1035 
1036  /* === Select a row index to use as the new pivot row =============== */
1037 
1038  pivot_row_length = pfree - pivot_row_start;
1039  if (pivot_row_length > 0) {
1040  /* pick the "pivot" row arbitrarily (first row in col) */
1041  pivot_row = A[Col[pivot_col].start];
1042  COLAMD_DEBUG3(("Pivotal row is %d\n", pivot_row));
1043  } else {
1044  /* there is no pivot row, since it is of zero length */
1045  pivot_row = Empty;
1046  COLAMD_ASSERT(pivot_row_length == 0);
1047  }
1048  COLAMD_ASSERT(Col[pivot_col].length > 0 || pivot_row_length == 0);
1049 
1050  /* === Approximate degree computation =============================== */
1051 
1052  /* Here begins the computation of the approximate degree. The column */
1053  /* score is the sum of the pivot row "length", plus the size of the */
1054  /* set differences of each row in the column minus the pattern of the */
1055  /* pivot row itself. The column ("thickness") itself is also */
1056  /* excluded from the column score (we thus use an approximate */
1057  /* external degree). */
1058 
1059  /* The time taken by the following code (compute set differences, and */
1060  /* add them up) is proportional to the size of the data structure */
1061  /* being scanned - that is, the sum of the sizes of each column in */
1062  /* the pivot row. Thus, the amortized time to compute a column score */
1063  /* is proportional to the size of that column (where size, in this */
1064  /* context, is the column "length", or the number of row indices */
1065  /* in that column). The number of row indices in a column is */
1066  /* monotonically non-decreasing, from the length of the original */
1067  /* column on input to colamd. */
1068 
1069  /* === Compute set differences ====================================== */
1070 
1071  COLAMD_DEBUG3(("** Computing set differences phase. **\n"));
1072 
1073  /* pivot row is currently dead - it will be revived later. */
1074 
1075  COLAMD_DEBUG3(("Pivot row: "));
1076  /* for each column in pivot row */
1077  rp = &A[pivot_row_start];
1078  rp_end = rp + pivot_row_length;
1079  while (rp < rp_end) {
1080  col = *rp++;
1081  COLAMD_ASSERT(Col[col].is_alive() && col != pivot_col);
1082  COLAMD_DEBUG3(("Col: %d\n", col));
1083 
1084  /* clear tags used to construct pivot row pattern */
1085  col_thickness = -Col[col].shared1.thickness;
1086  COLAMD_ASSERT(col_thickness > 0);
1087  Col[col].shared1.thickness = col_thickness;
1088 
1089  /* === Remove column from degree list =========================== */
1090 
1091  cur_score = Col[col].shared2.score;
1092  prev_col = Col[col].shared3.prev;
1093  next_col = Col[col].shared4.degree_next;
1094  COLAMD_ASSERT(cur_score >= 0);
1095  COLAMD_ASSERT(cur_score <= n_col);
1096  COLAMD_ASSERT(cur_score >= Empty);
1097  if (prev_col == Empty) {
1098  head[cur_score] = next_col;
1099  } else {
1100  Col[prev_col].shared4.degree_next = next_col;
1101  }
1102  if (next_col != Empty) {
1103  Col[next_col].shared3.prev = prev_col;
1104  }
1105 
1106  /* === Scan the column ========================================== */
1107 
1108  cp = &A[Col[col].start];
1109  cp_end = cp + Col[col].length;
1110  while (cp < cp_end) {
1111  /* get a row */
1112  row = *cp++;
1113  /* skip if dead */
1114  if (Row[row].is_dead()) {
1115  continue;
1116  }
1117  row_mark = Row[row].shared2.mark;
1118  COLAMD_ASSERT(row != pivot_row);
1119  set_difference = row_mark - tag_mark;
1120  /* check if the row has been seen yet */
1121  if (set_difference < 0) {
1122  COLAMD_ASSERT(Row[row].shared1.degree <= max_deg);
1123  set_difference = Row[row].shared1.degree;
1124  }
1125  /* subtract column thickness from this row's set difference */
1126  set_difference -= col_thickness;
1127  COLAMD_ASSERT(set_difference >= 0);
1128  /* absorb this row if the set difference becomes zero */
1129  if (set_difference == 0) {
1130  COLAMD_DEBUG3(("aggressive absorption. Row: %d\n", row));
1131  Row[row].kill();
1132  } else {
1133  /* save the new mark */
1134  Row[row].shared2.mark = set_difference + tag_mark;
1135  }
1136  }
1137  }
1138 
1139  /* === Add up set differences for each column ======================= */
1140 
1141  COLAMD_DEBUG3(("** Adding set differences phase. **\n"));
1142 
1143  /* for each column in pivot row */
1144  rp = &A[pivot_row_start];
1145  rp_end = rp + pivot_row_length;
1146  while (rp < rp_end) {
1147  /* get a column */
1148  col = *rp++;
1149  COLAMD_ASSERT(Col[col].is_alive() && col != pivot_col);
1150  hash = 0;
1151  cur_score = 0;
1152  cp = &A[Col[col].start];
1153  /* compact the column */
1154  new_cp = cp;
1155  cp_end = cp + Col[col].length;
1156 
1157  COLAMD_DEBUG4(("Adding set diffs for Col: %d.\n", col));
1158 
1159  while (cp < cp_end) {
1160  /* get a row */
1161  row = *cp++;
1162  COLAMD_ASSERT(row >= 0 && row < n_row);
1163  /* skip if dead */
1164  if (Row[row].is_dead()) {
1165  continue;
1166  }
1167  row_mark = Row[row].shared2.mark;
1168  COLAMD_ASSERT(row_mark > tag_mark);
1169  /* compact the column */
1170  *new_cp++ = row;
1171  /* compute hash function */
1172  hash += row;
1173  /* add set difference */
1174  cur_score += row_mark - tag_mark;
1175  /* integer overflow... */
1176  cur_score = numext::mini(cur_score, n_col);
1177  }
1178 
1179  /* recompute the column's length */
1180  Col[col].length = (IndexType)(new_cp - &A[Col[col].start]);
1181 
1182  /* === Further mass elimination ================================= */
1183 
1184  if (Col[col].length == 0) {
1185  COLAMD_DEBUG4(("further mass elimination. Col: %d\n", col));
1186  /* nothing left but the pivot row in this column */
1187  Col[col].kill_principal();
1188  pivot_row_degree -= Col[col].shared1.thickness;
1189  COLAMD_ASSERT(pivot_row_degree >= 0);
1190  /* order it */
1191  Col[col].shared2.order = k;
1192  /* increment order count by column thickness */
1193  k += Col[col].shared1.thickness;
1194  } else {
1195  /* === Prepare for supercolumn detection ==================== */
1196 
1197  COLAMD_DEBUG4(("Preparing supercol detection for Col: %d.\n", col));
1198 
1199  /* save score so far */
1200  Col[col].shared2.score = cur_score;
1201 
1202  /* add column to hash table, for supercolumn detection */
1203  hash %= n_col + 1;
1204 
1205  COLAMD_DEBUG4((" Hash = %d, n_col = %d.\n", hash, n_col));
1206  COLAMD_ASSERT(hash <= n_col);
1207 
1208  head_column = head[hash];
1209  if (head_column > Empty) {
1210  /* degree list "hash" is non-empty, use prev (shared3) of */
1211  /* first column in degree list as head of hash bucket */
1212  first_col = Col[head_column].shared3.headhash;
1213  Col[head_column].shared3.headhash = col;
1214  } else {
1215  /* degree list "hash" is empty, use head as hash bucket */
1216  first_col = -(head_column + 2);
1217  head[hash] = -(col + 2);
1218  }
1219  Col[col].shared4.hash_next = first_col;
1220 
1221  /* save hash function in Col [col].shared3.hash */
1222  Col[col].shared3.hash = (IndexType)hash;
1223  COLAMD_ASSERT(Col[col].is_alive());
1224  }
1225  }
1226 
1227  /* The approximate external column degree is now computed. */
1228 
1229  /* === Supercolumn detection ======================================== */
1230 
1231  COLAMD_DEBUG3(("** Supercolumn detection phase. **\n"));
1232 
1233  Colamd::detect_super_cols(Col, A, head, pivot_row_start, pivot_row_length);
1234 
1235  /* === Kill the pivotal column ====================================== */
1236 
1237  Col[pivot_col].kill_principal();
1238 
1239  /* === Clear mark =================================================== */
1240 
1241  tag_mark += (max_deg + 1);
1242  if (tag_mark >= max_mark) {
1243  COLAMD_DEBUG2(("clearing tag_mark\n"));
1244  tag_mark = Colamd::clear_mark(n_row, Row);
1245  }
1246 
1247  /* === Finalize the new pivot row, and column scores ================ */
1248 
1249  COLAMD_DEBUG3(("** Finalize scores phase. **\n"));
1250 
1251  /* for each column in pivot row */
1252  rp = &A[pivot_row_start];
1253  /* compact the pivot row */
1254  new_rp = rp;
1255  rp_end = rp + pivot_row_length;
1256  while (rp < rp_end) {
1257  col = *rp++;
1258  /* skip dead columns */
1259  if (Col[col].is_dead()) {
1260  continue;
1261  }
1262  *new_rp++ = col;
1263  /* add new pivot row to column */
1264  A[Col[col].start + (Col[col].length++)] = pivot_row;
1265 
1266  /* retrieve score so far and add on pivot row's degree. */
1267  /* (we wait until here for this in case the pivot */
1268  /* row's degree was reduced due to mass elimination). */
1269  cur_score = Col[col].shared2.score + pivot_row_degree;
1270 
1271  /* calculate the max possible score as the number of */
1272  /* external columns minus the 'k' value minus the */
1273  /* columns thickness */
1274  max_score = n_col - k - Col[col].shared1.thickness;
1275 
1276  /* make the score the external degree of the union-of-rows */
1277  cur_score -= Col[col].shared1.thickness;
1278 
1279  /* make sure score is less or equal than the max score */
1280  cur_score = numext::mini(cur_score, max_score);
1281  COLAMD_ASSERT(cur_score >= 0);
1282 
1283  /* store updated score */
1284  Col[col].shared2.score = cur_score;
1285 
1286  /* === Place column back in degree list ========================= */
1287 
1288  COLAMD_ASSERT(min_score >= 0);
1289  COLAMD_ASSERT(min_score <= n_col);
1290  COLAMD_ASSERT(cur_score >= 0);
1291  COLAMD_ASSERT(cur_score <= n_col);
1292  COLAMD_ASSERT(head[cur_score] >= Empty);
1293  next_col = head[cur_score];
1294  Col[col].shared4.degree_next = next_col;
1295  Col[col].shared3.prev = Empty;
1296  if (next_col != Empty) {
1297  Col[next_col].shared3.prev = col;
1298  }
1299  head[cur_score] = col;
1300 
1301  /* see if this score is less than current min */
1302  min_score = numext::mini(min_score, cur_score);
1303  }
1304 
1305  /* === Resurrect the new pivot row ================================== */
1306 
1307  if (pivot_row_degree > 0) {
1308  /* update pivot row length to reflect any cols that were killed */
1309  /* during super-col detection and mass elimination */
1310  Row[pivot_row].start = pivot_row_start;
1311  Row[pivot_row].length = (IndexType)(new_rp - &A[pivot_row_start]);
1312  Row[pivot_row].shared1.degree = pivot_row_degree;
1313  Row[pivot_row].shared2.mark = 0;
1314  /* pivot row is no longer dead */
1315  }
1316  }
1317 
1318  /* === All principal columns have now been ordered ====================== */
1319 
1320  return (ngarbage);
1321 }
1322 
1323 /* ========================================================================== */
1324 /* === order_children ======================================================= */
1325 /* ========================================================================== */
1326 
1327 /*
1328  The find_ordering routine has ordered all of the principal columns (the
1329  representatives of the supercolumns). The non-principal columns have not
1330  yet been ordered. This routine orders those columns by walking up the
1331  parent tree (a column is a child of the column which absorbed it). The
1332  final permutation vector is then placed in p [0 ... n_col-1], with p [0]
1333  being the first column, and p [n_col-1] being the last. It doesn't look
1334  like it at first glance, but be assured that this routine takes time linear
1335  in the number of columns. Although not immediately obvious, the time
1336  taken by this routine is O (n_col), that is, linear in the number of
1337  columns. Not user-callable.
1338 */
1339 template <typename IndexType>
1340 static inline void order_children(
1341  /* === Parameters ======================================================= */
1342 
1343  IndexType n_col, /* number of columns of A */
1344  ColStructure<IndexType> Col[], /* of size n_col+1 */
1345  IndexType p[] /* p [0 ... n_col-1] is the column permutation*/
1346 ) {
1347  /* === Local variables ================================================== */
1348 
1349  IndexType i; /* loop counter for all columns */
1350  IndexType c; /* column index */
1351  IndexType parent; /* index of column's parent */
1352  IndexType order; /* column's order */
1353 
1354  /* === Order each non-principal column ================================== */
1355 
1356  for (i = 0; i < n_col; i++) {
1357  /* find an un-ordered non-principal column */
1358  COLAMD_ASSERT(col_is_dead(Col, i));
1359  if (!Col[i].is_dead_principal() && Col[i].shared2.order == Empty) {
1360  parent = i;
1361  /* once found, find its principal parent */
1362  do {
1363  parent = Col[parent].shared1.parent;
1364  } while (!Col[parent].is_dead_principal());
1365 
1366  /* now, order all un-ordered non-principal columns along path */
1367  /* to this parent. collapse tree at the same time */
1368  c = i;
1369  /* get order of parent */
1370  order = Col[parent].shared2.order;
1371 
1372  do {
1373  COLAMD_ASSERT(Col[c].shared2.order == Empty);
1374 
1375  /* order this column */
1376  Col[c].shared2.order = order++;
1377  /* collapse tree */
1378  Col[c].shared1.parent = parent;
1379 
1380  /* get immediate parent of this column */
1381  c = Col[c].shared1.parent;
1382 
1383  /* continue until we hit an ordered column. There are */
1384  /* guaranteed not to be anymore unordered columns */
1385  /* above an ordered column */
1386  } while (Col[c].shared2.order == Empty);
1387 
1388  /* re-order the super_col parent to largest order for this group */
1389  Col[parent].shared2.order = order;
1390  }
1391  }
1392 
1393  /* === Generate the permutation ========================================= */
1394 
1395  for (c = 0; c < n_col; c++) {
1396  p[Col[c].shared2.order] = c;
1397  }
1398 }
1399 
1400 /* ========================================================================== */
1401 /* === detect_super_cols ==================================================== */
1402 /* ========================================================================== */
1403 
1404 /*
1405  Detects supercolumns by finding matches between columns in the hash buckets.
1406  Check amongst columns in the set A [row_start ... row_start + row_length-1].
1407  The columns under consideration are currently *not* in the degree lists,
1408  and have already been placed in the hash buckets.
1409 
1410  The hash bucket for columns whose hash function is equal to h is stored
1411  as follows:
1412 
1413  if head [h] is >= 0, then head [h] contains a degree list, so:
1414 
1415  head [h] is the first column in degree bucket h.
1416  Col [head [h]].headhash gives the first column in hash bucket h.
1417 
1418  otherwise, the degree list is empty, and:
1419 
1420  -(head [h] + 2) is the first column in hash bucket h.
1421 
1422  For a column c in a hash bucket, Col [c].shared3.prev is NOT a "previous
1423  column" pointer. Col [c].shared3.hash is used instead as the hash number
1424  for that column. The value of Col [c].shared4.hash_next is the next column
1425  in the same hash bucket.
1426 
1427  Assuming no, or "few" hash collisions, the time taken by this routine is
1428  linear in the sum of the sizes (lengths) of each column whose score has
1429  just been computed in the approximate degree computation.
1430  Not user-callable.
1431 */
1432 template <typename IndexType>
1433 static void detect_super_cols(
1434  /* === Parameters ======================================================= */
1435 
1436  ColStructure<IndexType> Col[], /* of size n_col+1 */
1437  IndexType A[], /* row indices of A */
1438  IndexType head[], /* head of degree lists and hash buckets */
1439  IndexType row_start, /* pointer to set of columns to check */
1440  IndexType row_length /* number of columns to check */
1441 ) {
1442  /* === Local variables ================================================== */
1443 
1444  IndexType hash; /* hash value for a column */
1445  IndexType *rp; /* pointer to a row */
1446  IndexType c; /* a column index */
1447  IndexType super_c; /* column index of the column to absorb into */
1448  IndexType *cp1; /* column pointer for column super_c */
1449  IndexType *cp2; /* column pointer for column c */
1450  IndexType length; /* length of column super_c */
1451  IndexType prev_c; /* column preceding c in hash bucket */
1452  IndexType i; /* loop counter */
1453  IndexType *rp_end; /* pointer to the end of the row */
1454  IndexType col; /* a column index in the row to check */
1455  IndexType head_column; /* first column in hash bucket or degree list */
1456  IndexType first_col; /* first column in hash bucket */
1457 
1458  /* === Consider each column in the row ================================== */
1459 
1460  rp = &A[row_start];
1461  rp_end = rp + row_length;
1462  while (rp < rp_end) {
1463  col = *rp++;
1464  if (Col[col].is_dead()) {
1465  continue;
1466  }
1467 
1468  /* get hash number for this column */
1469  hash = Col[col].shared3.hash;
1470  COLAMD_ASSERT(hash <= n_col);
1471 
1472  /* === Get the first column in this hash bucket ===================== */
1473 
1474  head_column = head[hash];
1475  if (head_column > Empty) {
1476  first_col = Col[head_column].shared3.headhash;
1477  } else {
1478  first_col = -(head_column + 2);
1479  }
1480 
1481  /* === Consider each column in the hash bucket ====================== */
1482 
1483  for (super_c = first_col; super_c != Empty; super_c = Col[super_c].shared4.hash_next) {
1484  COLAMD_ASSERT(Col[super_c].is_alive());
1485  COLAMD_ASSERT(Col[super_c].shared3.hash == hash);
1486  length = Col[super_c].length;
1487 
1488  /* prev_c is the column preceding column c in the hash bucket */
1489  prev_c = super_c;
1490 
1491  /* === Compare super_c with all columns after it ================ */
1492 
1493  for (c = Col[super_c].shared4.hash_next; c != Empty; c = Col[c].shared4.hash_next) {
1494  COLAMD_ASSERT(c != super_c);
1495  COLAMD_ASSERT(Col[c].is_alive());
1496  COLAMD_ASSERT(Col[c].shared3.hash == hash);
1497 
1498  /* not identical if lengths or scores are different */
1499  if (Col[c].length != length || Col[c].shared2.score != Col[super_c].shared2.score) {
1500  prev_c = c;
1501  continue;
1502  }
1503 
1504  /* compare the two columns */
1505  cp1 = &A[Col[super_c].start];
1506  cp2 = &A[Col[c].start];
1507 
1508  for (i = 0; i < length; i++) {
1509  /* the columns are "clean" (no dead rows) */
1510  COLAMD_ASSERT(cp1->is_alive());
1511  COLAMD_ASSERT(cp2->is_alive());
1512  /* row indices will same order for both supercols, */
1513  /* no gather scatter necessary */
1514  if (*cp1++ != *cp2++) {
1515  break;
1516  }
1517  }
1518 
1519  /* the two columns are different if the for-loop "broke" */
1520  if (i != length) {
1521  prev_c = c;
1522  continue;
1523  }
1524 
1525  /* === Got it! two columns are identical =================== */
1526 
1527  COLAMD_ASSERT(Col[c].shared2.score == Col[super_c].shared2.score);
1528 
1529  Col[super_c].shared1.thickness += Col[c].shared1.thickness;
1530  Col[c].shared1.parent = super_c;
1531  Col[c].kill_non_principal();
1532  /* order c later, in order_children() */
1533  Col[c].shared2.order = Empty;
1534  /* remove c from hash bucket */
1535  Col[prev_c].shared4.hash_next = Col[c].shared4.hash_next;
1536  }
1537  }
1538 
1539  /* === Empty this hash bucket ======================================= */
1540 
1541  if (head_column > Empty) {
1542  /* corresponding degree list "hash" is not empty */
1543  Col[head_column].shared3.headhash = Empty;
1544  } else {
1545  /* corresponding degree list "hash" is empty */
1546  head[hash] = Empty;
1547  }
1548  }
1549 }
1550 
1551 /* ========================================================================== */
1552 /* === garbage_collection =================================================== */
1553 /* ========================================================================== */
1554 
1555 /*
1556  Defragments and compacts columns and rows in the workspace A. Used when
1557  all available memory has been used while performing row merging. Returns
1558  the index of the first free position in A, after garbage collection. The
1559  time taken by this routine is linear is the size of the array A, which is
1560  itself linear in the number of nonzeros in the input matrix.
1561  Not user-callable.
1562 */
1563 template <typename IndexType>
1564 static IndexType garbage_collection /* returns the new value of pfree */
1565  (
1566  /* === Parameters ======================================================= */
1567 
1568  IndexType n_row, /* number of rows */
1569  IndexType n_col, /* number of columns */
1570  RowStructure<IndexType> Row[], /* row info */
1571  ColStructure<IndexType> Col[], /* column info */
1572  IndexType A[], /* A [0 ... Alen-1] holds the matrix */
1573  IndexType *pfree /* &A [0] ... pfree is in use */
1574  ) {
1575  /* === Local variables ================================================== */
1576 
1577  IndexType *psrc; /* source pointer */
1578  IndexType *pdest; /* destination pointer */
1579  IndexType j; /* counter */
1580  IndexType r; /* a row index */
1581  IndexType c; /* a column index */
1582  IndexType length; /* length of a row or column */
1583 
1584  /* === Defragment the columns =========================================== */
1585 
1586  pdest = &A[0];
1587  for (c = 0; c < n_col; c++) {
1588  if (Col[c].is_alive()) {
1589  psrc = &A[Col[c].start];
1590 
1591  /* move and compact the column */
1592  COLAMD_ASSERT(pdest <= psrc);
1593  Col[c].start = (IndexType)(pdest - &A[0]);
1594  length = Col[c].length;
1595  for (j = 0; j < length; j++) {
1596  r = *psrc++;
1597  if (Row[r].is_alive()) {
1598  *pdest++ = r;
1599  }
1600  }
1601  Col[c].length = (IndexType)(pdest - &A[Col[c].start]);
1602  }
1603  }
1604 
1605  /* === Prepare to defragment the rows =================================== */
1606 
1607  for (r = 0; r < n_row; r++) {
1608  if (Row[r].is_alive()) {
1609  if (Row[r].length == 0) {
1610  /* this row is of zero length. cannot compact it, so kill it */
1611  COLAMD_DEBUG3(("Defrag row kill\n"));
1612  Row[r].kill();
1613  } else {
1614  /* save first column index in Row [r].shared2.first_column */
1615  psrc = &A[Row[r].start];
1616  Row[r].shared2.first_column = *psrc;
1617  COLAMD_ASSERT(Row[r].is_alive());
1618  /* flag the start of the row with the one's complement of row */
1619  *psrc = ones_complement(r);
1620  }
1621  }
1622  }
1623 
1624  /* === Defragment the rows ============================================== */
1625 
1626  psrc = pdest;
1627  while (psrc < pfree) {
1628  /* find a negative number ... the start of a row */
1629  if (*psrc++ < 0) {
1630  psrc--;
1631  /* get the row index */
1632  r = ones_complement(*psrc);
1633  COLAMD_ASSERT(r >= 0 && r < n_row);
1634  /* restore first column index */
1635  *psrc = Row[r].shared2.first_column;
1636  COLAMD_ASSERT(Row[r].is_alive());
1637 
1638  /* move and compact the row */
1639  COLAMD_ASSERT(pdest <= psrc);
1640  Row[r].start = (IndexType)(pdest - &A[0]);
1641  length = Row[r].length;
1642  for (j = 0; j < length; j++) {
1643  c = *psrc++;
1644  if (Col[c].is_alive()) {
1645  *pdest++ = c;
1646  }
1647  }
1648  Row[r].length = (IndexType)(pdest - &A[Row[r].start]);
1649  }
1650  }
1651  /* ensure we found all the rows */
1652  COLAMD_ASSERT(debug_rows == 0);
1653 
1654  /* === Return the new value of pfree ==================================== */
1655 
1656  return ((IndexType)(pdest - &A[0]));
1657 }
1658 
1659 /* ========================================================================== */
1660 /* === clear_mark =========================================================== */
1661 /* ========================================================================== */
1662 
1663 /*
1664  Clears the Row [].shared2.mark array, and returns the new tag_mark.
1665  Return value is the new tag_mark. Not user-callable.
1666 */
1667 template <typename IndexType>
1668 static inline IndexType clear_mark /* return the new value for tag_mark */
1669  (
1670  /* === Parameters ======================================================= */
1671 
1672  IndexType n_row, /* number of rows in A */
1673  RowStructure<IndexType> Row[] /* Row [0 ... n_row-1].shared2.mark is set to zero */
1674  ) {
1675  /* === Local variables ================================================== */
1676 
1677  IndexType r;
1678 
1679  for (r = 0; r < n_row; r++) {
1680  if (Row[r].is_alive()) {
1681  Row[r].shared2.mark = 0;
1682  }
1683  }
1684  return (1);
1685 }
1686 
1687 } // namespace Colamd
1688 
1689 } // namespace internal
1690 #endif
int i
Definition: BiCGSTAB_step_by_step.cpp:9
#define COLAMD_ASSERT(expression)
Definition: Eigen_Colamd.h:250
#define COLAMD_DEBUG0(params)
Definition: Eigen_Colamd.h:244
#define COLAMD_DEBUG2(params)
Definition: Eigen_Colamd.h:246
#define COLAMD_DEBUG3(params)
Definition: Eigen_Colamd.h:247
#define COLAMD_DEBUG1(params)
Definition: Eigen_Colamd.h:245
#define COLAMD_DEBUG4(params)
Definition: Eigen_Colamd.h:248
m col(1)
m row(1)
float * p
Definition: Tutorial_Map_using.cpp:9
The matrix class, also used for vectors and row-vectors.
Definition: Eigen/Eigen/src/Core/Matrix.h:186
Data class to send an empty class over MPI.
Definition: MpiDataClass.h:110
char char char int int * k
Definition: level2_impl.h:374
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T maxi(const T &x, const T &y)
Definition: MathFunctions.h:926
EIGEN_DEVICE_FUNC EIGEN_ALWAYS_INLINE T mini(const T &x, const T &y)
Definition: MathFunctions.h:920
r
Definition: UniformPSDSelfTest.py:20
int c
Definition: calibrate.py:100
static IndexType garbage_collection(IndexType n_row, IndexType n_col, RowStructure< IndexType > Row[], ColStructure< IndexType > Col[], IndexType A[], IndexType *pfree)
Definition: Eigen_Colamd.h:1565
static void init_scoring(IndexType n_row, IndexType n_col, RowStructure< IndexType > Row[], ColStructure< IndexType > Col[], IndexType A[], IndexType head[], double knobs[NKnobs], IndexType *p_n_row2, IndexType *p_n_col2, IndexType *p_max_deg)
Definition: Eigen_Colamd.h:668
static IndexType find_ordering(IndexType n_row, IndexType n_col, IndexType Alen, RowStructure< IndexType > Row[], ColStructure< IndexType > Col[], IndexType A[], IndexType head[], IndexType n_col2, IndexType max_deg, IndexType pfree)
Definition: Eigen_Colamd.h:875
IndexType ones_complement(const IndexType r)
Definition: Eigen_Colamd.h:109
const int NKnobs
Definition: Eigen_Colamd.h:63
Status
Definition: Eigen_Colamd.h:89
@ ErrorRowIndexOutOfBounds
Definition: Eigen_Colamd.h:100
@ ErrorColLengthNegative
Definition: Eigen_Colamd.h:99
@ ErrorInternalError
Definition: Eigen_Colamd.h:102
@ ErrorANotPresent
Definition: Eigen_Colamd.h:92
@ ErrorNnzNegative
Definition: Eigen_Colamd.h:96
@ ErrorPNotPresent
Definition: Eigen_Colamd.h:93
@ ErrorNcolNegative
Definition: Eigen_Colamd.h:95
@ ErrorOutOfMemory
Definition: Eigen_Colamd.h:101
@ ErrorP0Nonzero
Definition: Eigen_Colamd.h:97
@ Ok
Definition: Eigen_Colamd.h:90
@ OkButJumbled
Definition: Eigen_Colamd.h:91
@ ErrorNrowNegative
Definition: Eigen_Colamd.h:94
@ ErrorATooSmall
Definition: Eigen_Colamd.h:98
static IndexType init_rows_cols(IndexType n_row, IndexType n_col, RowStructure< IndexType > Row[], ColStructure< IndexType > col[], IndexType A[], IndexType p[], IndexType stats[NStats])
Definition: Eigen_Colamd.h:479
static void order_children(IndexType n_col, ColStructure< IndexType > Col[], IndexType p[])
Definition: Eigen_Colamd.h:1340
IndexType colamd_r(IndexType n_row)
Definition: Eigen_Colamd.h:209
IndexType recommended(IndexType nnz, IndexType n_row, IndexType n_col)
Returns the recommended value of Alen.
Definition: Eigen_Colamd.h:267
static IndexType clear_mark(IndexType n_row, RowStructure< IndexType > Row[])
Definition: Eigen_Colamd.h:1669
KnobsStatsIndex
Definition: Eigen_Colamd.h:69
@ DenseRow
Definition: Eigen_Colamd.h:71
@ DefragCount
Definition: Eigen_Colamd.h:77
@ DenseCol
Definition: Eigen_Colamd.h:74
@ Info1
Definition: Eigen_Colamd.h:83
@ Info3
Definition: Eigen_Colamd.h:85
@ Info2
Definition: Eigen_Colamd.h:84
ColumnStatus
Definition: Eigen_Colamd.h:120
@ DeadNonPrincipal
Definition: Eigen_Colamd.h:120
@ DeadPrincipal
Definition: Eigen_Colamd.h:120
const int Empty
Definition: Eigen_Colamd.h:114
IndexType colamd_c(IndexType n_col)
Definition: Eigen_Colamd.h:204
static bool compute_ordering(IndexType n_row, IndexType n_col, IndexType Alen, IndexType *A, IndexType *p, double knobs[NKnobs], IndexType stats[NStats])
Computes a column ordering using the column approximate minimum degree ordering.
Definition: Eigen_Colamd.h:328
const int NStats
Definition: Eigen_Colamd.h:66
RowColumnStatus
Definition: Eigen_Colamd.h:117
@ Dead
Definition: Eigen_Colamd.h:117
@ Alive
Definition: Eigen_Colamd.h:117
static void set_defaults(double knobs[NKnobs])
set default parameters The use of this routine is optional.
Definition: Eigen_Colamd.h:295
static void detect_super_cols(ColStructure< IndexType > Col[], IndexType A[], IndexType head[], IndexType row_start, IndexType row_length)
Definition: Eigen_Colamd.h:1433
Definition: Eigen_Colamd.h:49
void start(const unsigned &i)
(Re-)start i-th timer
Definition: oomph_utilities.cc:243
Definition: Eigen_Colamd.h:128
union internal::Colamd::ColStructure::@779 shared1
bool is_alive() const
Definition: Eigen_Colamd.h:156
bool is_dead_principal() const
Definition: Eigen_Colamd.h:158
IndexType score
Definition: Eigen_Colamd.h:139
IndexType order
Definition: Eigen_Colamd.h:140
IndexType prev
Definition: Eigen_Colamd.h:146
bool is_dead() const
Definition: Eigen_Colamd.h:154
IndexType start
Definition: Eigen_Colamd.h:129
union internal::Colamd::ColStructure::@780 shared2
void kill_non_principal()
Definition: Eigen_Colamd.h:162
IndexType parent
Definition: Eigen_Colamd.h:135
IndexType hash
Definition: Eigen_Colamd.h:145
IndexType thickness
Definition: Eigen_Colamd.h:133
IndexType hash_next
Definition: Eigen_Colamd.h:151
void kill_principal()
Definition: Eigen_Colamd.h:160
IndexType degree_next
Definition: Eigen_Colamd.h:150
union internal::Colamd::ColStructure::@782 shared4
union internal::Colamd::ColStructure::@781 shared3
IndexType headhash
Definition: Eigen_Colamd.h:143
IndexType length
Definition: Eigen_Colamd.h:131
Definition: Eigen_Colamd.h:166
bool is_alive() const
Definition: Eigen_Colamd.h:180
IndexType length
Definition: Eigen_Colamd.h:168
IndexType first_column
Definition: Eigen_Colamd.h:175
IndexType degree
Definition: Eigen_Colamd.h:170
bool is_dead() const
Definition: Eigen_Colamd.h:178
union internal::Colamd::RowStructure::@784 shared2
union internal::Colamd::RowStructure::@783 shared1
IndexType mark
Definition: Eigen_Colamd.h:174
void kill()
Definition: Eigen_Colamd.h:182
IndexType start
Definition: Eigen_Colamd.h:167
IndexType p
Definition: Eigen_Colamd.h:171
std::ptrdiff_t j
Definition: tut_arithmetic_redux_minmax.cpp:2