Back to home page

EIC code displayed by LXR

 
 

    


Warning, file /include/eigen3/Eigen/src/OrderingMethods/Eigen_Colamd.h was not indexed or was modified since last indexation (in which case cross-reference links may be missing, inaccurate or erroneous).

0001 // // This file is part of Eigen, a lightweight C++ template library
0002 // for linear algebra.
0003 //
0004 // Copyright (C) 2012 Desire Nuentsa Wakam <desire.nuentsa_wakam@inria.fr>
0005 //
0006 // This Source Code Form is subject to the terms of the Mozilla
0007 // Public License v. 2.0. If a copy of the MPL was not distributed
0008 // with this file, You can obtain one at http://mozilla.org/MPL/2.0/.
0009 
0010 // This file is modified from the colamd/symamd library. The copyright is below
0011 
0012 //   The authors of the code itself are Stefan I. Larimore and Timothy A.
0013 //   Davis (davis@cise.ufl.edu), University of Florida.  The algorithm was
0014 //   developed in collaboration with John Gilbert, Xerox PARC, and Esmond
0015 //   Ng, Oak Ridge National Laboratory.
0016 //
0017 //     Date:
0018 //
0019 //   September 8, 2003.  Version 2.3.
0020 //
0021 //     Acknowledgements:
0022 //
0023 //   This work was supported by the National Science Foundation, under
0024 //   grants DMS-9504974 and DMS-9803599.
0025 //
0026 //     Notice:
0027 //
0028 //   Copyright (c) 1998-2003 by the University of Florida.
0029 //   All Rights Reserved.
0030 //
0031 //   THIS MATERIAL IS PROVIDED AS IS, WITH ABSOLUTELY NO WARRANTY
0032 //   EXPRESSED OR IMPLIED.  ANY USE IS AT YOUR OWN RISK.
0033 //
0034 //   Permission is hereby granted to use, copy, modify, and/or distribute
0035 //   this program, provided that the Copyright, this License, and the
0036 //   Availability of the original version is retained on all copies and made
0037 //   accessible to the end-user of any code or package that includes COLAMD
0038 //   or any modified version of COLAMD.
0039 //
0040 //     Availability:
0041 //
0042 //   The colamd/symamd library is available at
0043 //
0044 //       http://www.suitesparse.com
0045 
0046 
0047 #ifndef EIGEN_COLAMD_H
0048 #define EIGEN_COLAMD_H
0049 
0050 namespace internal {
0051 
0052 namespace Colamd {
0053 
0054 /* Ensure that debugging is turned off: */
0055 #ifndef COLAMD_NDEBUG
0056 #define COLAMD_NDEBUG
0057 #endif /* NDEBUG */
0058 
0059 
0060 /* ========================================================================== */
0061 /* === Knob and statistics definitions ====================================== */
0062 /* ========================================================================== */
0063 
0064 /* size of the knobs [ ] array.  Only knobs [0..1] are currently used. */
0065 const int NKnobs = 20;
0066 
0067 /* number of output statistics.  Only stats [0..6] are currently used. */
0068 const int NStats = 20;
0069 
0070 /* Indices into knobs and stats array. */
0071 enum KnobsStatsIndex {
0072   /* knobs [0] and stats [0]: dense row knob and output statistic. */
0073   DenseRow = 0,
0074 
0075   /* knobs [1] and stats [1]: dense column knob and output statistic. */
0076   DenseCol = 1,
0077 
0078   /* stats [2]: memory defragmentation count output statistic */
0079   DefragCount = 2,
0080 
0081   /* stats [3]: colamd status:  zero OK, > 0 warning or notice, < 0 error */
0082   Status = 3,
0083 
0084   /* stats [4..6]: error info, or info on jumbled columns */
0085   Info1 = 4,
0086   Info2 = 5,
0087   Info3 = 6
0088 };
0089 
0090 /* error codes returned in stats [3]: */
0091 enum Status {
0092   Ok = 0,
0093   OkButJumbled = 1,
0094   ErrorANotPresent = -1,
0095   ErrorPNotPresent = -2,
0096   ErrorNrowNegative = -3,
0097   ErrorNcolNegative = -4,
0098   ErrorNnzNegative = -5,
0099   ErrorP0Nonzero = -6,
0100   ErrorATooSmall = -7,
0101   ErrorColLengthNegative = -8,
0102   ErrorRowIndexOutOfBounds = -9,
0103   ErrorOutOfMemory = -10,
0104   ErrorInternalError = -999
0105 };
0106 /* ========================================================================== */
0107 /* === Definitions ========================================================== */
0108 /* ========================================================================== */
0109 
0110 template <typename IndexType>
0111 IndexType ones_complement(const IndexType r) {
0112   return (-(r)-1);
0113 }
0114 
0115 /* -------------------------------------------------------------------------- */
0116 const int Empty = -1;
0117 
0118 /* Row and column status */
0119 enum RowColumnStatus {
0120   Alive = 0,
0121   Dead = -1
0122 };
0123 
0124 /* Column status */
0125 enum ColumnStatus {
0126   DeadPrincipal = -1,
0127   DeadNonPrincipal = -2
0128 };
0129 
0130 /* ========================================================================== */
0131 /* === Colamd reporting mechanism =========================================== */
0132 /* ========================================================================== */
0133 
0134 // == Row and Column structures ==
0135 template <typename IndexType>
0136 struct ColStructure
0137 {
0138   IndexType start ;   /* index for A of first row in this column, or Dead */
0139   /* if column is dead */
0140   IndexType length ;  /* number of rows in this column */
0141   union
0142   {
0143     IndexType thickness ; /* number of original columns represented by this */
0144     /* col, if the column is alive */
0145     IndexType parent ;  /* parent in parent tree super-column structure, if */
0146     /* the column is dead */
0147   } shared1 ;
0148   union
0149   {
0150     IndexType score ; /* the score used to maintain heap, if col is alive */
0151     IndexType order ; /* pivot ordering of this column, if col is dead */
0152   } shared2 ;
0153   union
0154   {
0155     IndexType headhash ;  /* head of a hash bucket, if col is at the head of */
0156     /* a degree list */
0157     IndexType hash ;  /* hash value, if col is not in a degree list */
0158     IndexType prev ;  /* previous column in degree list, if col is in a */
0159     /* degree list (but not at the head of a degree list) */
0160   } shared3 ;
0161   union
0162   {
0163     IndexType degree_next ; /* next column, if col is in a degree list */
0164     IndexType hash_next ;   /* next column, if col is in a hash list */
0165   } shared4 ;
0166 
0167   inline bool is_dead() const { return start < Alive; }
0168 
0169   inline bool is_alive() const { return start >= Alive; }
0170 
0171   inline bool is_dead_principal() const { return start == DeadPrincipal; }
0172 
0173   inline void kill_principal() { start = DeadPrincipal; }
0174 
0175   inline void kill_non_principal() { start = DeadNonPrincipal; }
0176 
0177 };
0178 
0179 template <typename IndexType>
0180 struct RowStructure
0181 {
0182   IndexType start ;   /* index for A of first col in this row */
0183   IndexType length ;  /* number of principal columns in this row */
0184   union
0185   {
0186     IndexType degree ;  /* number of principal & non-principal columns in row */
0187     IndexType p ;   /* used as a row pointer in init_rows_cols () */
0188   } shared1 ;
0189   union
0190   {
0191     IndexType mark ;  /* for computing set differences and marking dead rows*/
0192     IndexType first_column ;/* first column in row (used in garbage collection) */
0193   } shared2 ;
0194 
0195   inline bool is_dead() const { return shared2.mark < Alive; }
0196 
0197   inline bool is_alive() const { return shared2.mark >= Alive; }
0198 
0199   inline void kill() { shared2.mark = Dead; }
0200 
0201 };
0202 
0203 /* ========================================================================== */
0204 /* === Colamd recommended memory size ======================================= */
0205 /* ========================================================================== */
0206 
0207 /*
0208   The recommended length Alen of the array A passed to colamd is given by
0209   the COLAMD_RECOMMENDED (nnz, n_row, n_col) macro.  It returns -1 if any
0210   argument is negative.  2*nnz space is required for the row and column
0211   indices of the matrix. colamd_c (n_col) + colamd_r (n_row) space is
0212   required for the Col and Row arrays, respectively, which are internal to
0213   colamd.  An additional n_col space is the minimal amount of "elbow room",
0214   and nnz/5 more space is recommended for run time efficiency.
0215 
0216   This macro is not needed when using symamd.
0217 
0218   Explicit typecast to IndexType added Sept. 23, 2002, COLAMD version 2.2, to avoid
0219   gcc -pedantic warning messages.
0220 */
0221 template <typename IndexType>
0222 inline IndexType colamd_c(IndexType n_col)
0223 { return IndexType( ((n_col) + 1) * sizeof (ColStructure<IndexType>) / sizeof (IndexType) ) ; }
0224 
0225 template <typename IndexType>
0226 inline IndexType  colamd_r(IndexType n_row)
0227 { return IndexType(((n_row) + 1) * sizeof (RowStructure<IndexType>) / sizeof (IndexType)); }
0228 
0229 // Prototypes of non-user callable routines
0230 template <typename IndexType>
0231 static IndexType init_rows_cols (IndexType n_row, IndexType n_col, RowStructure<IndexType> Row [], ColStructure<IndexType> col [], IndexType A [], IndexType p [], IndexType stats[NStats] );
0232 
0233 template <typename IndexType>
0234 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);
0235 
0236 template <typename IndexType>
0237 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);
0238 
0239 template <typename IndexType>
0240 static void order_children (IndexType n_col, ColStructure<IndexType> Col [], IndexType p []);
0241 
0242 template <typename IndexType>
0243 static void detect_super_cols (ColStructure<IndexType> Col [], IndexType A [], IndexType head [], IndexType row_start, IndexType row_length ) ;
0244 
0245 template <typename IndexType>
0246 static IndexType garbage_collection (IndexType n_row, IndexType n_col, RowStructure<IndexType> Row [], ColStructure<IndexType> Col [], IndexType A [], IndexType *pfree) ;
0247 
0248 template <typename IndexType>
0249 static inline  IndexType clear_mark (IndexType n_row, RowStructure<IndexType> Row [] ) ;
0250 
0251 /* === No debugging ========================================================= */
0252 
0253 #define COLAMD_DEBUG0(params) ;
0254 #define COLAMD_DEBUG1(params) ;
0255 #define COLAMD_DEBUG2(params) ;
0256 #define COLAMD_DEBUG3(params) ;
0257 #define COLAMD_DEBUG4(params) ;
0258 
0259 #define COLAMD_ASSERT(expression) ((void) 0)
0260 
0261 
0262 /**
0263  * \brief Returns the recommended value of Alen
0264  *
0265  * Returns recommended value of Alen for use by colamd.
0266  * Returns -1 if any input argument is negative.
0267  * The use of this routine or macro is optional.
0268  * Note that the macro uses its arguments   more than once,
0269  * so be careful for side effects, if you pass expressions as arguments to COLAMD_RECOMMENDED.
0270  *
0271  * \param nnz nonzeros in A
0272  * \param n_row number of rows in A
0273  * \param n_col number of columns in A
0274  * \return recommended value of Alen for use by colamd
0275  */
0276 template <typename IndexType>
0277 inline IndexType recommended ( IndexType nnz, IndexType n_row, IndexType n_col)
0278 {
0279   if ((nnz) < 0 || (n_row) < 0 || (n_col) < 0)
0280     return (-1);
0281   else
0282     return (2 * (nnz) + colamd_c (n_col) + colamd_r (n_row) + (n_col) + ((nnz) / 5));
0283 }
0284 
0285 /**
0286  * \brief set default parameters  The use of this routine is optional.
0287  *
0288  * Colamd: rows with more than (knobs [DenseRow] * n_col)
0289  * entries are removed prior to ordering.  Columns with more than
0290  * (knobs [DenseCol] * n_row) entries are removed prior to
0291  * ordering, and placed last in the output column ordering.
0292  *
0293  * DenseRow and DenseCol are defined as 0 and 1,
0294  * respectively, in colamd.h.  Default values of these two knobs
0295  * are both 0.5.  Currently, only knobs [0] and knobs [1] are
0296  * used, but future versions may use more knobs.  If so, they will
0297  * be properly set to their defaults by the future version of
0298  * colamd_set_defaults, so that the code that calls colamd will
0299  * not need to change, assuming that you either use
0300  * colamd_set_defaults, or pass a (double *) NULL pointer as the
0301  * knobs array to colamd or symamd.
0302  *
0303  * \param knobs parameter settings for colamd
0304  */
0305 
0306 static inline void set_defaults(double knobs[NKnobs])
0307 {
0308   /* === Local variables ================================================== */
0309 
0310   int i ;
0311 
0312   if (!knobs)
0313   {
0314     return ;      /* no knobs to initialize */
0315   }
0316   for (i = 0 ; i < NKnobs ; i++)
0317   {
0318     knobs [i] = 0 ;
0319   }
0320   knobs [Colamd::DenseRow] = 0.5 ;  /* ignore rows over 50% dense */
0321   knobs [Colamd::DenseCol] = 0.5 ;  /* ignore columns over 50% dense */
0322 }
0323 
0324 /**
0325  * \brief  Computes a column ordering using the column approximate minimum degree ordering
0326  *
0327  * Computes a column ordering (Q) of A such that P(AQ)=LU or
0328  * (AQ)'AQ=LL' have less fill-in and require fewer floating point
0329  * operations than factorizing the unpermuted matrix A or A'A,
0330  * respectively.
0331  *
0332  *
0333  * \param n_row number of rows in A
0334  * \param n_col number of columns in A
0335  * \param Alen, size of the array A
0336  * \param A row indices of the matrix, of size ALen
0337  * \param p column pointers of A, of size n_col+1
0338  * \param knobs parameter settings for colamd
0339  * \param stats colamd output statistics and error codes
0340  */
0341 template <typename IndexType>
0342 static bool compute_ordering(IndexType n_row, IndexType n_col, IndexType Alen, IndexType *A, IndexType *p, double knobs[NKnobs], IndexType stats[NStats])
0343 {
0344   /* === Local variables ================================================== */
0345 
0346   IndexType i ;     /* loop index */
0347   IndexType nnz ;     /* nonzeros in A */
0348   IndexType Row_size ;    /* size of Row [], in integers */
0349   IndexType Col_size ;    /* size of Col [], in integers */
0350   IndexType need ;      /* minimum required length of A */
0351   Colamd::RowStructure<IndexType> *Row ;   /* pointer into A of Row [0..n_row] array */
0352   Colamd::ColStructure<IndexType> *Col ;   /* pointer into A of Col [0..n_col] array */
0353   IndexType n_col2 ;    /* number of non-dense, non-empty columns */
0354   IndexType n_row2 ;    /* number of non-dense, non-empty rows */
0355   IndexType ngarbage ;    /* number of garbage collections performed */
0356   IndexType max_deg ;   /* maximum row degree */
0357   double default_knobs [NKnobs] ; /* default knobs array */
0358 
0359 
0360   /* === Check the input arguments ======================================== */
0361 
0362   if (!stats)
0363   {
0364     COLAMD_DEBUG0 (("colamd: stats not present\n")) ;
0365     return (false) ;
0366   }
0367   for (i = 0 ; i < NStats ; i++)
0368   {
0369     stats [i] = 0 ;
0370   }
0371   stats [Colamd::Status] = Colamd::Ok ;
0372   stats [Colamd::Info1] = -1 ;
0373   stats [Colamd::Info2] = -1 ;
0374 
0375   if (!A)   /* A is not present */
0376   {
0377     stats [Colamd::Status] = Colamd::ErrorANotPresent ;
0378     COLAMD_DEBUG0 (("colamd: A not present\n")) ;
0379     return (false) ;
0380   }
0381 
0382   if (!p)   /* p is not present */
0383   {
0384     stats [Colamd::Status] = Colamd::ErrorPNotPresent ;
0385     COLAMD_DEBUG0 (("colamd: p not present\n")) ;
0386     return (false) ;
0387   }
0388 
0389   if (n_row < 0)  /* n_row must be >= 0 */
0390   {
0391     stats [Colamd::Status] = Colamd::ErrorNrowNegative ;
0392     stats [Colamd::Info1] = n_row ;
0393     COLAMD_DEBUG0 (("colamd: nrow negative %d\n", n_row)) ;
0394     return (false) ;
0395   }
0396 
0397   if (n_col < 0)  /* n_col must be >= 0 */
0398   {
0399     stats [Colamd::Status] = Colamd::ErrorNcolNegative ;
0400     stats [Colamd::Info1] = n_col ;
0401     COLAMD_DEBUG0 (("colamd: ncol negative %d\n", n_col)) ;
0402     return (false) ;
0403   }
0404 
0405   nnz = p [n_col] ;
0406   if (nnz < 0)  /* nnz must be >= 0 */
0407   {
0408     stats [Colamd::Status] = Colamd::ErrorNnzNegative ;
0409     stats [Colamd::Info1] = nnz ;
0410     COLAMD_DEBUG0 (("colamd: number of entries negative %d\n", nnz)) ;
0411     return (false) ;
0412   }
0413 
0414   if (p [0] != 0)
0415   {
0416     stats [Colamd::Status] = Colamd::ErrorP0Nonzero ;
0417     stats [Colamd::Info1] = p [0] ;
0418     COLAMD_DEBUG0 (("colamd: p[0] not zero %d\n", p [0])) ;
0419     return (false) ;
0420   }
0421 
0422   /* === If no knobs, set default knobs =================================== */
0423 
0424   if (!knobs)
0425   {
0426     set_defaults (default_knobs) ;
0427     knobs = default_knobs ;
0428   }
0429 
0430   /* === Allocate the Row and Col arrays from array A ===================== */
0431 
0432   Col_size = colamd_c (n_col) ;
0433   Row_size = colamd_r (n_row) ;
0434   need = 2*nnz + n_col + Col_size + Row_size ;
0435 
0436   if (need > Alen)
0437   {
0438     /* not enough space in array A to perform the ordering */
0439     stats [Colamd::Status] = Colamd::ErrorATooSmall ;
0440     stats [Colamd::Info1] = need ;
0441     stats [Colamd::Info2] = Alen ;
0442     COLAMD_DEBUG0 (("colamd: Need Alen >= %d, given only Alen = %d\n", need,Alen));
0443     return (false) ;
0444   }
0445 
0446   Alen -= Col_size + Row_size ;
0447   Col = (ColStructure<IndexType> *) &A [Alen] ;
0448   Row = (RowStructure<IndexType> *) &A [Alen + Col_size] ;
0449 
0450   /* === Construct the row and column data structures ===================== */
0451 
0452   if (!Colamd::init_rows_cols (n_row, n_col, Row, Col, A, p, stats))
0453   {
0454     /* input matrix is invalid */
0455     COLAMD_DEBUG0 (("colamd: Matrix invalid\n")) ;
0456     return (false) ;
0457   }
0458 
0459   /* === Initialize scores, kill dense rows/columns ======================= */
0460 
0461   Colamd::init_scoring (n_row, n_col, Row, Col, A, p, knobs,
0462         &n_row2, &n_col2, &max_deg) ;
0463 
0464   /* === Order the supercolumns =========================================== */
0465 
0466   ngarbage = Colamd::find_ordering (n_row, n_col, Alen, Row, Col, A, p,
0467                 n_col2, max_deg, 2*nnz) ;
0468 
0469   /* === Order the non-principal columns ================================== */
0470 
0471   Colamd::order_children (n_col, Col, p) ;
0472 
0473   /* === Return statistics in stats ======================================= */
0474 
0475   stats [Colamd::DenseRow] = n_row - n_row2 ;
0476   stats [Colamd::DenseCol] = n_col - n_col2 ;
0477   stats [Colamd::DefragCount] = ngarbage ;
0478   COLAMD_DEBUG0 (("colamd: done.\n")) ;
0479   return (true) ;
0480 }
0481 
0482 /* ========================================================================== */
0483 /* === NON-USER-CALLABLE ROUTINES: ========================================== */
0484 /* ========================================================================== */
0485 
0486 /* There are no user-callable routines beyond this point in the file */
0487 
0488 /* ========================================================================== */
0489 /* === init_rows_cols ======================================================= */
0490 /* ========================================================================== */
0491 
0492 /*
0493   Takes the column form of the matrix in A and creates the row form of the
0494   matrix.  Also, row and column attributes are stored in the Col and Row
0495   structs.  If the columns are un-sorted or contain duplicate row indices,
0496   this routine will also sort and remove duplicate row indices from the
0497   column form of the matrix.  Returns false if the matrix is invalid,
0498   true otherwise.  Not user-callable.
0499 */
0500 template <typename IndexType>
0501 static IndexType init_rows_cols  /* returns true if OK, or false otherwise */
0502   (
0503     /* === Parameters ======================================================= */
0504 
0505     IndexType n_row,      /* number of rows of A */
0506     IndexType n_col,      /* number of columns of A */
0507     RowStructure<IndexType> Row [],    /* of size n_row+1 */
0508     ColStructure<IndexType> Col [],    /* of size n_col+1 */
0509     IndexType A [],     /* row indices of A, of size Alen */
0510     IndexType p [],     /* pointers to columns in A, of size n_col+1 */
0511     IndexType stats [NStats]  /* colamd statistics */
0512     )
0513 {
0514   /* === Local variables ================================================== */
0515 
0516   IndexType col ;     /* a column index */
0517   IndexType row ;     /* a row index */
0518   IndexType *cp ;     /* a column pointer */
0519   IndexType *cp_end ;   /* a pointer to the end of a column */
0520   IndexType *rp ;     /* a row pointer */
0521   IndexType *rp_end ;   /* a pointer to the end of a row */
0522   IndexType last_row ;    /* previous row */
0523 
0524   /* === Initialize columns, and check column pointers ==================== */
0525 
0526   for (col = 0 ; col < n_col ; col++)
0527   {
0528     Col [col].start = p [col] ;
0529     Col [col].length = p [col+1] - p [col] ;
0530 
0531     if ((Col [col].length) < 0) // extra parentheses to work-around gcc bug 10200
0532     {
0533       /* column pointers must be non-decreasing */
0534       stats [Colamd::Status] = Colamd::ErrorColLengthNegative ;
0535       stats [Colamd::Info1] = col ;
0536       stats [Colamd::Info2] = Col [col].length ;
0537       COLAMD_DEBUG0 (("colamd: col %d length %d < 0\n", col, Col [col].length)) ;
0538       return (false) ;
0539     }
0540 
0541     Col [col].shared1.thickness = 1 ;
0542     Col [col].shared2.score = 0 ;
0543     Col [col].shared3.prev = Empty ;
0544     Col [col].shared4.degree_next = Empty ;
0545   }
0546 
0547   /* p [0..n_col] no longer needed, used as "head" in subsequent routines */
0548 
0549   /* === Scan columns, compute row degrees, and check row indices ========= */
0550 
0551   stats [Info3] = 0 ;  /* number of duplicate or unsorted row indices*/
0552 
0553   for (row = 0 ; row < n_row ; row++)
0554   {
0555     Row [row].length = 0 ;
0556     Row [row].shared2.mark = -1 ;
0557   }
0558 
0559   for (col = 0 ; col < n_col ; col++)
0560   {
0561     last_row = -1 ;
0562 
0563     cp = &A [p [col]] ;
0564     cp_end = &A [p [col+1]] ;
0565 
0566     while (cp < cp_end)
0567     {
0568       row = *cp++ ;
0569 
0570       /* make sure row indices within range */
0571       if (row < 0 || row >= n_row)
0572       {
0573     stats [Colamd::Status] = Colamd::ErrorRowIndexOutOfBounds ;
0574     stats [Colamd::Info1] = col ;
0575     stats [Colamd::Info2] = row ;
0576     stats [Colamd::Info3] = n_row ;
0577     COLAMD_DEBUG0 (("colamd: row %d col %d out of bounds\n", row, col)) ;
0578     return (false) ;
0579       }
0580 
0581       if (row <= last_row || Row [row].shared2.mark == col)
0582       {
0583     /* row index are unsorted or repeated (or both), thus col */
0584     /* is jumbled.  This is a notice, not an error condition. */
0585     stats [Colamd::Status] = Colamd::OkButJumbled ;
0586     stats [Colamd::Info1] = col ;
0587     stats [Colamd::Info2] = row ;
0588     (stats [Colamd::Info3]) ++ ;
0589     COLAMD_DEBUG1 (("colamd: row %d col %d unsorted/duplicate\n",row,col));
0590       }
0591 
0592       if (Row [row].shared2.mark != col)
0593       {
0594     Row [row].length++ ;
0595       }
0596       else
0597       {
0598     /* this is a repeated entry in the column, */
0599     /* it will be removed */
0600     Col [col].length-- ;
0601       }
0602 
0603       /* mark the row as having been seen in this column */
0604       Row [row].shared2.mark = col ;
0605 
0606       last_row = row ;
0607     }
0608   }
0609 
0610   /* === Compute row pointers ============================================= */
0611 
0612   /* row form of the matrix starts directly after the column */
0613   /* form of matrix in A */
0614   Row [0].start = p [n_col] ;
0615   Row [0].shared1.p = Row [0].start ;
0616   Row [0].shared2.mark = -1 ;
0617   for (row = 1 ; row < n_row ; row++)
0618   {
0619     Row [row].start = Row [row-1].start + Row [row-1].length ;
0620     Row [row].shared1.p = Row [row].start ;
0621     Row [row].shared2.mark = -1 ;
0622   }
0623 
0624   /* === Create row form ================================================== */
0625 
0626   if (stats [Status] == OkButJumbled)
0627   {
0628     /* if cols jumbled, watch for repeated row indices */
0629     for (col = 0 ; col < n_col ; col++)
0630     {
0631       cp = &A [p [col]] ;
0632       cp_end = &A [p [col+1]] ;
0633       while (cp < cp_end)
0634       {
0635     row = *cp++ ;
0636     if (Row [row].shared2.mark != col)
0637     {
0638       A [(Row [row].shared1.p)++] = col ;
0639       Row [row].shared2.mark = col ;
0640     }
0641       }
0642     }
0643   }
0644   else
0645   {
0646     /* if cols not jumbled, we don't need the mark (this is faster) */
0647     for (col = 0 ; col < n_col ; col++)
0648     {
0649       cp = &A [p [col]] ;
0650       cp_end = &A [p [col+1]] ;
0651       while (cp < cp_end)
0652       {
0653     A [(Row [*cp++].shared1.p)++] = col ;
0654       }
0655     }
0656   }
0657 
0658   /* === Clear the row marks and set row degrees ========================== */
0659 
0660   for (row = 0 ; row < n_row ; row++)
0661   {
0662     Row [row].shared2.mark = 0 ;
0663     Row [row].shared1.degree = Row [row].length ;
0664   }
0665 
0666   /* === See if we need to re-create columns ============================== */
0667 
0668   if (stats [Status] == OkButJumbled)
0669   {
0670     COLAMD_DEBUG0 (("colamd: reconstructing column form, matrix jumbled\n")) ;
0671 
0672 
0673     /* === Compute col pointers ========================================= */
0674 
0675     /* col form of the matrix starts at A [0]. */
0676     /* Note, we may have a gap between the col form and the row */
0677     /* form if there were duplicate entries, if so, it will be */
0678     /* removed upon the first garbage collection */
0679     Col [0].start = 0 ;
0680     p [0] = Col [0].start ;
0681     for (col = 1 ; col < n_col ; col++)
0682     {
0683       /* note that the lengths here are for pruned columns, i.e. */
0684       /* no duplicate row indices will exist for these columns */
0685       Col [col].start = Col [col-1].start + Col [col-1].length ;
0686       p [col] = Col [col].start ;
0687     }
0688 
0689     /* === Re-create col form =========================================== */
0690 
0691     for (row = 0 ; row < n_row ; row++)
0692     {
0693       rp = &A [Row [row].start] ;
0694       rp_end = rp + Row [row].length ;
0695       while (rp < rp_end)
0696       {
0697     A [(p [*rp++])++] = row ;
0698       }
0699     }
0700   }
0701 
0702   /* === Done.  Matrix is not (or no longer) jumbled ====================== */
0703 
0704   return (true) ;
0705 }
0706 
0707 
0708 /* ========================================================================== */
0709 /* === init_scoring ========================================================= */
0710 /* ========================================================================== */
0711 
0712 /*
0713   Kills dense or empty columns and rows, calculates an initial score for
0714   each column, and places all columns in the degree lists.  Not user-callable.
0715 */
0716 template <typename IndexType>
0717 static void init_scoring
0718   (
0719     /* === Parameters ======================================================= */
0720 
0721     IndexType n_row,      /* number of rows of A */
0722     IndexType n_col,      /* number of columns of A */
0723     RowStructure<IndexType> Row [],    /* of size n_row+1 */
0724     ColStructure<IndexType> Col [],    /* of size n_col+1 */
0725     IndexType A [],     /* column form and row form of A */
0726     IndexType head [],    /* of size n_col+1 */
0727     double knobs [NKnobs],/* parameters */
0728     IndexType *p_n_row2,    /* number of non-dense, non-empty rows */
0729     IndexType *p_n_col2,    /* number of non-dense, non-empty columns */
0730     IndexType *p_max_deg    /* maximum row degree */
0731     )
0732 {
0733   /* === Local variables ================================================== */
0734 
0735   IndexType c ;     /* a column index */
0736   IndexType r, row ;    /* a row index */
0737   IndexType *cp ;     /* a column pointer */
0738   IndexType deg ;     /* degree of a row or column */
0739   IndexType *cp_end ;   /* a pointer to the end of a column */
0740   IndexType *new_cp ;   /* new column pointer */
0741   IndexType col_length ;    /* length of pruned column */
0742   IndexType score ;     /* current column score */
0743   IndexType n_col2 ;    /* number of non-dense, non-empty columns */
0744   IndexType n_row2 ;    /* number of non-dense, non-empty rows */
0745   IndexType dense_row_count ; /* remove rows with more entries than this */
0746   IndexType dense_col_count ; /* remove cols with more entries than this */
0747   IndexType min_score ;   /* smallest column score */
0748   IndexType max_deg ;   /* maximum row degree */
0749   IndexType next_col ;    /* Used to add to degree list.*/
0750 
0751 
0752   /* === Extract knobs ==================================================== */
0753 
0754   dense_row_count = numext::maxi(IndexType(0), numext::mini(IndexType(knobs [Colamd::DenseRow] * n_col), n_col)) ;
0755   dense_col_count = numext::maxi(IndexType(0), numext::mini(IndexType(knobs [Colamd::DenseCol] * n_row), n_row)) ;
0756   COLAMD_DEBUG1 (("colamd: densecount: %d %d\n", dense_row_count, dense_col_count)) ;
0757   max_deg = 0 ;
0758   n_col2 = n_col ;
0759   n_row2 = n_row ;
0760 
0761   /* === Kill empty columns =============================================== */
0762 
0763   /* Put the empty columns at the end in their natural order, so that LU */
0764   /* factorization can proceed as far as possible. */
0765   for (c = n_col-1 ; c >= 0 ; c--)
0766   {
0767     deg = Col [c].length ;
0768     if (deg == 0)
0769     {
0770       /* this is a empty column, kill and order it last */
0771       Col [c].shared2.order = --n_col2 ;
0772       Col[c].kill_principal() ;
0773     }
0774   }
0775   COLAMD_DEBUG1 (("colamd: null columns killed: %d\n", n_col - n_col2)) ;
0776 
0777   /* === Kill dense columns =============================================== */
0778 
0779   /* Put the dense columns at the end, in their natural order */
0780   for (c = n_col-1 ; c >= 0 ; c--)
0781   {
0782     /* skip any dead columns */
0783     if (Col[c].is_dead())
0784     {
0785       continue ;
0786     }
0787     deg = Col [c].length ;
0788     if (deg > dense_col_count)
0789     {
0790       /* this is a dense column, kill and order it last */
0791       Col [c].shared2.order = --n_col2 ;
0792       /* decrement the row degrees */
0793       cp = &A [Col [c].start] ;
0794       cp_end = cp + Col [c].length ;
0795       while (cp < cp_end)
0796       {
0797     Row [*cp++].shared1.degree-- ;
0798       }
0799       Col[c].kill_principal() ;
0800     }
0801   }
0802   COLAMD_DEBUG1 (("colamd: Dense and null columns killed: %d\n", n_col - n_col2)) ;
0803 
0804   /* === Kill dense and empty rows ======================================== */
0805 
0806   for (r = 0 ; r < n_row ; r++)
0807   {
0808     deg = Row [r].shared1.degree ;
0809     COLAMD_ASSERT (deg >= 0 && deg <= n_col) ;
0810     if (deg > dense_row_count || deg == 0)
0811     {
0812       /* kill a dense or empty row */
0813       Row[r].kill() ;
0814       --n_row2 ;
0815     }
0816     else
0817     {
0818       /* keep track of max degree of remaining rows */
0819       max_deg = numext::maxi(max_deg, deg) ;
0820     }
0821   }
0822   COLAMD_DEBUG1 (("colamd: Dense and null rows killed: %d\n", n_row - n_row2)) ;
0823 
0824   /* === Compute initial column scores ==================================== */
0825 
0826   /* At this point the row degrees are accurate.  They reflect the number */
0827   /* of "live" (non-dense) columns in each row.  No empty rows exist. */
0828   /* Some "live" columns may contain only dead rows, however.  These are */
0829   /* pruned in the code below. */
0830 
0831   /* now find the initial matlab score for each column */
0832   for (c = n_col-1 ; c >= 0 ; c--)
0833   {
0834     /* skip dead column */
0835     if (Col[c].is_dead())
0836     {
0837       continue ;
0838     }
0839     score = 0 ;
0840     cp = &A [Col [c].start] ;
0841     new_cp = cp ;
0842     cp_end = cp + Col [c].length ;
0843     while (cp < cp_end)
0844     {
0845       /* get a row */
0846       row = *cp++ ;
0847       /* skip if dead */
0848       if (Row[row].is_dead())
0849       {
0850     continue ;
0851       }
0852       /* compact the column */
0853       *new_cp++ = row ;
0854       /* add row's external degree */
0855       score += Row [row].shared1.degree - 1 ;
0856       /* guard against integer overflow */
0857       score = numext::mini(score, n_col) ;
0858     }
0859     /* determine pruned column length */
0860     col_length = (IndexType) (new_cp - &A [Col [c].start]) ;
0861     if (col_length == 0)
0862     {
0863       /* a newly-made null column (all rows in this col are "dense" */
0864       /* and have already been killed) */
0865       COLAMD_DEBUG2 (("Newly null killed: %d\n", c)) ;
0866       Col [c].shared2.order = --n_col2 ;
0867       Col[c].kill_principal() ;
0868     }
0869     else
0870     {
0871       /* set column length and set score */
0872       COLAMD_ASSERT (score >= 0) ;
0873       COLAMD_ASSERT (score <= n_col) ;
0874       Col [c].length = col_length ;
0875       Col [c].shared2.score = score ;
0876     }
0877   }
0878   COLAMD_DEBUG1 (("colamd: Dense, null, and newly-null columns killed: %d\n",
0879           n_col-n_col2)) ;
0880 
0881   /* At this point, all empty rows and columns are dead.  All live columns */
0882   /* are "clean" (containing no dead rows) and simplicial (no supercolumns */
0883   /* yet).  Rows may contain dead columns, but all live rows contain at */
0884   /* least one live column. */
0885 
0886   /* === Initialize degree lists ========================================== */
0887 
0888 
0889   /* clear the hash buckets */
0890   for (c = 0 ; c <= n_col ; c++)
0891   {
0892     head [c] = Empty ;
0893   }
0894   min_score = n_col ;
0895   /* place in reverse order, so low column indices are at the front */
0896   /* of the lists.  This is to encourage natural tie-breaking */
0897   for (c = n_col-1 ; c >= 0 ; c--)
0898   {
0899     /* only add principal columns to degree lists */
0900     if (Col[c].is_alive())
0901     {
0902       COLAMD_DEBUG4 (("place %d score %d minscore %d ncol %d\n",
0903               c, Col [c].shared2.score, min_score, n_col)) ;
0904 
0905       /* === Add columns score to DList =============================== */
0906 
0907       score = Col [c].shared2.score ;
0908 
0909       COLAMD_ASSERT (min_score >= 0) ;
0910       COLAMD_ASSERT (min_score <= n_col) ;
0911       COLAMD_ASSERT (score >= 0) ;
0912       COLAMD_ASSERT (score <= n_col) ;
0913       COLAMD_ASSERT (head [score] >= Empty) ;
0914 
0915       /* now add this column to dList at proper score location */
0916       next_col = head [score] ;
0917       Col [c].shared3.prev = Empty ;
0918       Col [c].shared4.degree_next = next_col ;
0919 
0920       /* if there already was a column with the same score, set its */
0921       /* previous pointer to this new column */
0922       if (next_col != Empty)
0923       {
0924     Col [next_col].shared3.prev = c ;
0925       }
0926       head [score] = c ;
0927 
0928       /* see if this score is less than current min */
0929       min_score = numext::mini(min_score, score) ;
0930 
0931 
0932     }
0933   }
0934 
0935 
0936   /* === Return number of remaining columns, and max row degree =========== */
0937 
0938   *p_n_col2 = n_col2 ;
0939   *p_n_row2 = n_row2 ;
0940   *p_max_deg = max_deg ;
0941 }
0942 
0943 
0944 /* ========================================================================== */
0945 /* === find_ordering ======================================================== */
0946 /* ========================================================================== */
0947 
0948 /*
0949   Order the principal columns of the supercolumn form of the matrix
0950   (no supercolumns on input).  Uses a minimum approximate column minimum
0951   degree ordering method.  Not user-callable.
0952 */
0953 template <typename IndexType>
0954 static IndexType find_ordering /* return the number of garbage collections */
0955   (
0956     /* === Parameters ======================================================= */
0957 
0958     IndexType n_row,      /* number of rows of A */
0959     IndexType n_col,      /* number of columns of A */
0960     IndexType Alen,     /* size of A, 2*nnz + n_col or larger */
0961     RowStructure<IndexType> Row [],    /* of size n_row+1 */
0962     ColStructure<IndexType> Col [],    /* of size n_col+1 */
0963     IndexType A [],     /* column form and row form of A */
0964     IndexType head [],    /* of size n_col+1 */
0965     IndexType n_col2,     /* Remaining columns to order */
0966     IndexType max_deg,    /* Maximum row degree */
0967     IndexType pfree     /* index of first free slot (2*nnz on entry) */
0968     )
0969 {
0970   /* === Local variables ================================================== */
0971 
0972   IndexType k ;     /* current pivot ordering step */
0973   IndexType pivot_col ;   /* current pivot column */
0974   IndexType *cp ;     /* a column pointer */
0975   IndexType *rp ;     /* a row pointer */
0976   IndexType pivot_row ;   /* current pivot row */
0977   IndexType *new_cp ;   /* modified column pointer */
0978   IndexType *new_rp ;   /* modified row pointer */
0979   IndexType pivot_row_start ; /* pointer to start of pivot row */
0980   IndexType pivot_row_degree ;  /* number of columns in pivot row */
0981   IndexType pivot_row_length ;  /* number of supercolumns in pivot row */
0982   IndexType pivot_col_score ; /* score of pivot column */
0983   IndexType needed_memory ;   /* free space needed for pivot row */
0984   IndexType *cp_end ;   /* pointer to the end of a column */
0985   IndexType *rp_end ;   /* pointer to the end of a row */
0986   IndexType row ;     /* a row index */
0987   IndexType col ;     /* a column index */
0988   IndexType max_score ;   /* maximum possible score */
0989   IndexType cur_score ;   /* score of current column */
0990   unsigned int hash ;   /* hash value for supernode detection */
0991   IndexType head_column ;   /* head of hash bucket */
0992   IndexType first_col ;   /* first column in hash bucket */
0993   IndexType tag_mark ;    /* marker value for mark array */
0994   IndexType row_mark ;    /* Row [row].shared2.mark */
0995   IndexType set_difference ;  /* set difference size of row with pivot row */
0996   IndexType min_score ;   /* smallest column score */
0997   IndexType col_thickness ;   /* "thickness" (no. of columns in a supercol) */
0998   IndexType max_mark ;    /* maximum value of tag_mark */
0999   IndexType pivot_col_thickness ; /* number of columns represented by pivot col */
1000   IndexType prev_col ;    /* Used by Dlist operations. */
1001   IndexType next_col ;    /* Used by Dlist operations. */
1002   IndexType ngarbage ;    /* number of garbage collections performed */
1003 
1004 
1005   /* === Initialization and clear mark ==================================== */
1006 
1007   max_mark = INT_MAX - n_col ;  /* INT_MAX defined in <limits.h> */
1008   tag_mark = Colamd::clear_mark (n_row, Row) ;
1009   min_score = 0 ;
1010   ngarbage = 0 ;
1011   COLAMD_DEBUG1 (("colamd: Ordering, n_col2=%d\n", n_col2)) ;
1012 
1013   /* === Order the columns ================================================ */
1014 
1015   for (k = 0 ; k < n_col2 ; /* 'k' is incremented below */)
1016   {
1017 
1018     /* === Select pivot column, and order it ============================ */
1019 
1020     /* make sure degree list isn't empty */
1021     COLAMD_ASSERT (min_score >= 0) ;
1022     COLAMD_ASSERT (min_score <= n_col) ;
1023     COLAMD_ASSERT (head [min_score] >= Empty) ;
1024 
1025     /* get pivot column from head of minimum degree list */
1026     while (min_score < n_col && head [min_score] == Empty)
1027     {
1028       min_score++ ;
1029     }
1030     pivot_col = head [min_score] ;
1031     COLAMD_ASSERT (pivot_col >= 0 && pivot_col <= n_col) ;
1032     next_col = Col [pivot_col].shared4.degree_next ;
1033     head [min_score] = next_col ;
1034     if (next_col != Empty)
1035     {
1036       Col [next_col].shared3.prev = Empty ;
1037     }
1038 
1039     COLAMD_ASSERT (Col[pivot_col].is_alive()) ;
1040     COLAMD_DEBUG3 (("Pivot col: %d\n", pivot_col)) ;
1041 
1042     /* remember score for defrag check */
1043     pivot_col_score = Col [pivot_col].shared2.score ;
1044 
1045     /* the pivot column is the kth column in the pivot order */
1046     Col [pivot_col].shared2.order = k ;
1047 
1048     /* increment order count by column thickness */
1049     pivot_col_thickness = Col [pivot_col].shared1.thickness ;
1050     k += pivot_col_thickness ;
1051     COLAMD_ASSERT (pivot_col_thickness > 0) ;
1052 
1053     /* === Garbage_collection, if necessary ============================= */
1054 
1055     needed_memory = numext::mini(pivot_col_score, n_col - k) ;
1056     if (pfree + needed_memory >= Alen)
1057     {
1058       pfree = Colamd::garbage_collection (n_row, n_col, Row, Col, A, &A [pfree]) ;
1059       ngarbage++ ;
1060       /* after garbage collection we will have enough */
1061       COLAMD_ASSERT (pfree + needed_memory < Alen) ;
1062       /* garbage collection has wiped out the Row[].shared2.mark array */
1063       tag_mark = Colamd::clear_mark (n_row, Row) ;
1064 
1065     }
1066 
1067     /* === Compute pivot row pattern ==================================== */
1068 
1069     /* get starting location for this new merged row */
1070     pivot_row_start = pfree ;
1071 
1072     /* initialize new row counts to zero */
1073     pivot_row_degree = 0 ;
1074 
1075     /* tag pivot column as having been visited so it isn't included */
1076     /* in merged pivot row */
1077     Col [pivot_col].shared1.thickness = -pivot_col_thickness ;
1078 
1079     /* pivot row is the union of all rows in the pivot column pattern */
1080     cp = &A [Col [pivot_col].start] ;
1081     cp_end = cp + Col [pivot_col].length ;
1082     while (cp < cp_end)
1083     {
1084       /* get a row */
1085       row = *cp++ ;
1086       COLAMD_DEBUG4 (("Pivot col pattern %d %d\n", Row[row].is_alive(), row)) ;
1087       /* skip if row is dead */
1088       if (Row[row].is_dead())
1089       {
1090     continue ;
1091       }
1092       rp = &A [Row [row].start] ;
1093       rp_end = rp + Row [row].length ;
1094       while (rp < rp_end)
1095       {
1096     /* get a column */
1097     col = *rp++ ;
1098     /* add the column, if alive and untagged */
1099     col_thickness = Col [col].shared1.thickness ;
1100     if (col_thickness > 0 && Col[col].is_alive())
1101     {
1102       /* tag column in pivot row */
1103       Col [col].shared1.thickness = -col_thickness ;
1104       COLAMD_ASSERT (pfree < Alen) ;
1105       /* place column in pivot row */
1106       A [pfree++] = col ;
1107       pivot_row_degree += col_thickness ;
1108     }
1109       }
1110     }
1111 
1112     /* clear tag on pivot column */
1113     Col [pivot_col].shared1.thickness = pivot_col_thickness ;
1114     max_deg = numext::maxi(max_deg, pivot_row_degree) ;
1115 
1116 
1117     /* === Kill all rows used to construct pivot row ==================== */
1118 
1119     /* also kill pivot row, temporarily */
1120     cp = &A [Col [pivot_col].start] ;
1121     cp_end = cp + Col [pivot_col].length ;
1122     while (cp < cp_end)
1123     {
1124       /* may be killing an already dead row */
1125       row = *cp++ ;
1126       COLAMD_DEBUG3 (("Kill row in pivot col: %d\n", row)) ;
1127       Row[row].kill() ;
1128     }
1129 
1130     /* === Select a row index to use as the new pivot row =============== */
1131 
1132     pivot_row_length = pfree - pivot_row_start ;
1133     if (pivot_row_length > 0)
1134     {
1135       /* pick the "pivot" row arbitrarily (first row in col) */
1136       pivot_row = A [Col [pivot_col].start] ;
1137       COLAMD_DEBUG3 (("Pivotal row is %d\n", pivot_row)) ;
1138     }
1139     else
1140     {
1141       /* there is no pivot row, since it is of zero length */
1142       pivot_row = Empty ;
1143       COLAMD_ASSERT (pivot_row_length == 0) ;
1144     }
1145     COLAMD_ASSERT (Col [pivot_col].length > 0 || pivot_row_length == 0) ;
1146 
1147     /* === Approximate degree computation =============================== */
1148 
1149     /* Here begins the computation of the approximate degree.  The column */
1150     /* score is the sum of the pivot row "length", plus the size of the */
1151     /* set differences of each row in the column minus the pattern of the */
1152     /* pivot row itself.  The column ("thickness") itself is also */
1153     /* excluded from the column score (we thus use an approximate */
1154     /* external degree). */
1155 
1156     /* The time taken by the following code (compute set differences, and */
1157     /* add them up) is proportional to the size of the data structure */
1158     /* being scanned - that is, the sum of the sizes of each column in */
1159     /* the pivot row.  Thus, the amortized time to compute a column score */
1160     /* is proportional to the size of that column (where size, in this */
1161     /* context, is the column "length", or the number of row indices */
1162     /* in that column).  The number of row indices in a column is */
1163     /* monotonically non-decreasing, from the length of the original */
1164     /* column on input to colamd. */
1165 
1166     /* === Compute set differences ====================================== */
1167 
1168     COLAMD_DEBUG3 (("** Computing set differences phase. **\n")) ;
1169 
1170     /* pivot row is currently dead - it will be revived later. */
1171 
1172     COLAMD_DEBUG3 (("Pivot row: ")) ;
1173     /* for each column in pivot row */
1174     rp = &A [pivot_row_start] ;
1175     rp_end = rp + pivot_row_length ;
1176     while (rp < rp_end)
1177     {
1178       col = *rp++ ;
1179       COLAMD_ASSERT (Col[col].is_alive() && col != pivot_col) ;
1180       COLAMD_DEBUG3 (("Col: %d\n", col)) ;
1181 
1182       /* clear tags used to construct pivot row pattern */
1183       col_thickness = -Col [col].shared1.thickness ;
1184       COLAMD_ASSERT (col_thickness > 0) ;
1185       Col [col].shared1.thickness = col_thickness ;
1186 
1187       /* === Remove column from degree list =========================== */
1188 
1189       cur_score = Col [col].shared2.score ;
1190       prev_col = Col [col].shared3.prev ;
1191       next_col = Col [col].shared4.degree_next ;
1192       COLAMD_ASSERT (cur_score >= 0) ;
1193       COLAMD_ASSERT (cur_score <= n_col) ;
1194       COLAMD_ASSERT (cur_score >= Empty) ;
1195       if (prev_col == Empty)
1196       {
1197     head [cur_score] = next_col ;
1198       }
1199       else
1200       {
1201     Col [prev_col].shared4.degree_next = next_col ;
1202       }
1203       if (next_col != Empty)
1204       {
1205     Col [next_col].shared3.prev = prev_col ;
1206       }
1207 
1208       /* === Scan the column ========================================== */
1209 
1210       cp = &A [Col [col].start] ;
1211       cp_end = cp + Col [col].length ;
1212       while (cp < cp_end)
1213       {
1214     /* get a row */
1215     row = *cp++ ;
1216     /* skip if dead */
1217     if (Row[row].is_dead())
1218     {
1219       continue ;
1220     }
1221   row_mark = Row [row].shared2.mark ;
1222     COLAMD_ASSERT (row != pivot_row) ;
1223     set_difference = row_mark - tag_mark ;
1224     /* check if the row has been seen yet */
1225     if (set_difference < 0)
1226     {
1227       COLAMD_ASSERT (Row [row].shared1.degree <= max_deg) ;
1228       set_difference = Row [row].shared1.degree ;
1229     }
1230     /* subtract column thickness from this row's set difference */
1231     set_difference -= col_thickness ;
1232     COLAMD_ASSERT (set_difference >= 0) ;
1233     /* absorb this row if the set difference becomes zero */
1234     if (set_difference == 0)
1235     {
1236       COLAMD_DEBUG3 (("aggressive absorption. Row: %d\n", row)) ;
1237       Row[row].kill() ;
1238     }
1239     else
1240     {
1241       /* save the new mark */
1242       Row [row].shared2.mark = set_difference + tag_mark ;
1243     }
1244       }
1245     }
1246 
1247 
1248     /* === Add up set differences for each column ======================= */
1249 
1250     COLAMD_DEBUG3 (("** Adding set differences phase. **\n")) ;
1251 
1252     /* for each column in pivot row */
1253     rp = &A [pivot_row_start] ;
1254     rp_end = rp + pivot_row_length ;
1255     while (rp < rp_end)
1256     {
1257       /* get a column */
1258       col = *rp++ ;
1259       COLAMD_ASSERT (Col[col].is_alive() && col != pivot_col) ;
1260       hash = 0 ;
1261       cur_score = 0 ;
1262       cp = &A [Col [col].start] ;
1263       /* compact the column */
1264       new_cp = cp ;
1265       cp_end = cp + Col [col].length ;
1266 
1267       COLAMD_DEBUG4 (("Adding set diffs for Col: %d.\n", col)) ;
1268 
1269       while (cp < cp_end)
1270       {
1271     /* get a row */
1272     row = *cp++ ;
1273     COLAMD_ASSERT(row >= 0 && row < n_row) ;
1274     /* skip if dead */
1275     if (Row [row].is_dead())
1276     {
1277       continue ;
1278     }
1279   row_mark = Row [row].shared2.mark ;
1280     COLAMD_ASSERT (row_mark > tag_mark) ;
1281     /* compact the column */
1282     *new_cp++ = row ;
1283     /* compute hash function */
1284     hash += row ;
1285     /* add set difference */
1286     cur_score += row_mark - tag_mark ;
1287     /* integer overflow... */
1288     cur_score = numext::mini(cur_score, n_col) ;
1289       }
1290 
1291       /* recompute the column's length */
1292       Col [col].length = (IndexType) (new_cp - &A [Col [col].start]) ;
1293 
1294       /* === Further mass elimination ================================= */
1295 
1296       if (Col [col].length == 0)
1297       {
1298     COLAMD_DEBUG4 (("further mass elimination. Col: %d\n", col)) ;
1299     /* nothing left but the pivot row in this column */
1300     Col[col].kill_principal() ;
1301     pivot_row_degree -= Col [col].shared1.thickness ;
1302     COLAMD_ASSERT (pivot_row_degree >= 0) ;
1303     /* order it */
1304     Col [col].shared2.order = k ;
1305     /* increment order count by column thickness */
1306     k += Col [col].shared1.thickness ;
1307       }
1308       else
1309       {
1310     /* === Prepare for supercolumn detection ==================== */
1311 
1312     COLAMD_DEBUG4 (("Preparing supercol detection for Col: %d.\n", col)) ;
1313 
1314     /* save score so far */
1315     Col [col].shared2.score = cur_score ;
1316 
1317     /* add column to hash table, for supercolumn detection */
1318     hash %= n_col + 1 ;
1319 
1320     COLAMD_DEBUG4 ((" Hash = %d, n_col = %d.\n", hash, n_col)) ;
1321     COLAMD_ASSERT (hash <= n_col) ;
1322 
1323     head_column = head [hash] ;
1324     if (head_column > Empty)
1325     {
1326       /* degree list "hash" is non-empty, use prev (shared3) of */
1327       /* first column in degree list as head of hash bucket */
1328       first_col = Col [head_column].shared3.headhash ;
1329       Col [head_column].shared3.headhash = col ;
1330     }
1331     else
1332     {
1333       /* degree list "hash" is empty, use head as hash bucket */
1334       first_col = - (head_column + 2) ;
1335       head [hash] = - (col + 2) ;
1336     }
1337     Col [col].shared4.hash_next = first_col ;
1338 
1339     /* save hash function in Col [col].shared3.hash */
1340     Col [col].shared3.hash = (IndexType) hash ;
1341     COLAMD_ASSERT (Col[col].is_alive()) ;
1342       }
1343     }
1344 
1345     /* The approximate external column degree is now computed.  */
1346 
1347     /* === Supercolumn detection ======================================== */
1348 
1349     COLAMD_DEBUG3 (("** Supercolumn detection phase. **\n")) ;
1350 
1351     Colamd::detect_super_cols (Col, A, head, pivot_row_start, pivot_row_length) ;
1352 
1353     /* === Kill the pivotal column ====================================== */
1354 
1355     Col[pivot_col].kill_principal() ;
1356 
1357     /* === Clear mark =================================================== */
1358 
1359     tag_mark += (max_deg + 1) ;
1360     if (tag_mark >= max_mark)
1361     {
1362       COLAMD_DEBUG2 (("clearing tag_mark\n")) ;
1363       tag_mark = Colamd::clear_mark (n_row, Row) ;
1364     }
1365 
1366     /* === Finalize the new pivot row, and column scores ================ */
1367 
1368     COLAMD_DEBUG3 (("** Finalize scores phase. **\n")) ;
1369 
1370     /* for each column in pivot row */
1371     rp = &A [pivot_row_start] ;
1372     /* compact the pivot row */
1373     new_rp = rp ;
1374     rp_end = rp + pivot_row_length ;
1375     while (rp < rp_end)
1376     {
1377       col = *rp++ ;
1378       /* skip dead columns */
1379       if (Col[col].is_dead())
1380       {
1381     continue ;
1382       }
1383       *new_rp++ = col ;
1384       /* add new pivot row to column */
1385       A [Col [col].start + (Col [col].length++)] = pivot_row ;
1386 
1387       /* retrieve score so far and add on pivot row's degree. */
1388       /* (we wait until here for this in case the pivot */
1389       /* row's degree was reduced due to mass elimination). */
1390       cur_score = Col [col].shared2.score + pivot_row_degree ;
1391 
1392       /* calculate the max possible score as the number of */
1393       /* external columns minus the 'k' value minus the */
1394       /* columns thickness */
1395       max_score = n_col - k - Col [col].shared1.thickness ;
1396 
1397       /* make the score the external degree of the union-of-rows */
1398       cur_score -= Col [col].shared1.thickness ;
1399 
1400       /* make sure score is less or equal than the max score */
1401       cur_score = numext::mini(cur_score, max_score) ;
1402       COLAMD_ASSERT (cur_score >= 0) ;
1403 
1404       /* store updated score */
1405       Col [col].shared2.score = cur_score ;
1406 
1407       /* === Place column back in degree list ========================= */
1408 
1409       COLAMD_ASSERT (min_score >= 0) ;
1410       COLAMD_ASSERT (min_score <= n_col) ;
1411       COLAMD_ASSERT (cur_score >= 0) ;
1412       COLAMD_ASSERT (cur_score <= n_col) ;
1413       COLAMD_ASSERT (head [cur_score] >= Empty) ;
1414       next_col = head [cur_score] ;
1415       Col [col].shared4.degree_next = next_col ;
1416       Col [col].shared3.prev = Empty ;
1417       if (next_col != Empty)
1418       {
1419     Col [next_col].shared3.prev = col ;
1420       }
1421       head [cur_score] = col ;
1422 
1423       /* see if this score is less than current min */
1424       min_score = numext::mini(min_score, cur_score) ;
1425 
1426     }
1427 
1428     /* === Resurrect the new pivot row ================================== */
1429 
1430     if (pivot_row_degree > 0)
1431     {
1432       /* update pivot row length to reflect any cols that were killed */
1433       /* during super-col detection and mass elimination */
1434       Row [pivot_row].start  = pivot_row_start ;
1435       Row [pivot_row].length = (IndexType) (new_rp - &A[pivot_row_start]) ;
1436       Row [pivot_row].shared1.degree = pivot_row_degree ;
1437       Row [pivot_row].shared2.mark = 0 ;
1438       /* pivot row is no longer dead */
1439     }
1440   }
1441 
1442   /* === All principal columns have now been ordered ====================== */
1443 
1444   return (ngarbage) ;
1445 }
1446 
1447 
1448 /* ========================================================================== */
1449 /* === order_children ======================================================= */
1450 /* ========================================================================== */
1451 
1452 /*
1453   The find_ordering routine has ordered all of the principal columns (the
1454   representatives of the supercolumns).  The non-principal columns have not
1455   yet been ordered.  This routine orders those columns by walking up the
1456   parent tree (a column is a child of the column which absorbed it).  The
1457   final permutation vector is then placed in p [0 ... n_col-1], with p [0]
1458   being the first column, and p [n_col-1] being the last.  It doesn't look
1459   like it at first glance, but be assured that this routine takes time linear
1460   in the number of columns.  Although not immediately obvious, the time
1461   taken by this routine is O (n_col), that is, linear in the number of
1462   columns.  Not user-callable.
1463 */
1464 template <typename IndexType>
1465 static inline  void order_children
1466 (
1467   /* === Parameters ======================================================= */
1468 
1469   IndexType n_col,      /* number of columns of A */
1470   ColStructure<IndexType> Col [],    /* of size n_col+1 */
1471   IndexType p []      /* p [0 ... n_col-1] is the column permutation*/
1472   )
1473 {
1474   /* === Local variables ================================================== */
1475 
1476   IndexType i ;     /* loop counter for all columns */
1477   IndexType c ;     /* column index */
1478   IndexType parent ;    /* index of column's parent */
1479   IndexType order ;     /* column's order */
1480 
1481   /* === Order each non-principal column ================================== */
1482 
1483   for (i = 0 ; i < n_col ; i++)
1484   {
1485     /* find an un-ordered non-principal column */
1486     COLAMD_ASSERT (col_is_dead(Col, i)) ;
1487     if (!Col[i].is_dead_principal() && Col [i].shared2.order == Empty)
1488     {
1489       parent = i ;
1490       /* once found, find its principal parent */
1491       do
1492       {
1493     parent = Col [parent].shared1.parent ;
1494       } while (!Col[parent].is_dead_principal()) ;
1495 
1496       /* now, order all un-ordered non-principal columns along path */
1497       /* to this parent.  collapse tree at the same time */
1498       c = i ;
1499       /* get order of parent */
1500       order = Col [parent].shared2.order ;
1501 
1502       do
1503       {
1504     COLAMD_ASSERT (Col [c].shared2.order == Empty) ;
1505 
1506     /* order this column */
1507     Col [c].shared2.order = order++ ;
1508     /* collaps tree */
1509     Col [c].shared1.parent = parent ;
1510 
1511     /* get immediate parent of this column */
1512     c = Col [c].shared1.parent ;
1513 
1514     /* continue until we hit an ordered column.  There are */
1515     /* guaranteed not to be anymore unordered columns */
1516     /* above an ordered column */
1517       } while (Col [c].shared2.order == Empty) ;
1518 
1519       /* re-order the super_col parent to largest order for this group */
1520       Col [parent].shared2.order = order ;
1521     }
1522   }
1523 
1524   /* === Generate the permutation ========================================= */
1525 
1526   for (c = 0 ; c < n_col ; c++)
1527   {
1528     p [Col [c].shared2.order] = c ;
1529   }
1530 }
1531 
1532 
1533 /* ========================================================================== */
1534 /* === detect_super_cols ==================================================== */
1535 /* ========================================================================== */
1536 
1537 /*
1538   Detects supercolumns by finding matches between columns in the hash buckets.
1539   Check amongst columns in the set A [row_start ... row_start + row_length-1].
1540   The columns under consideration are currently *not* in the degree lists,
1541   and have already been placed in the hash buckets.
1542 
1543   The hash bucket for columns whose hash function is equal to h is stored
1544   as follows:
1545 
1546   if head [h] is >= 0, then head [h] contains a degree list, so:
1547 
1548   head [h] is the first column in degree bucket h.
1549   Col [head [h]].headhash gives the first column in hash bucket h.
1550 
1551   otherwise, the degree list is empty, and:
1552 
1553   -(head [h] + 2) is the first column in hash bucket h.
1554 
1555   For a column c in a hash bucket, Col [c].shared3.prev is NOT a "previous
1556   column" pointer.  Col [c].shared3.hash is used instead as the hash number
1557   for that column.  The value of Col [c].shared4.hash_next is the next column
1558   in the same hash bucket.
1559 
1560   Assuming no, or "few" hash collisions, the time taken by this routine is
1561   linear in the sum of the sizes (lengths) of each column whose score has
1562   just been computed in the approximate degree computation.
1563   Not user-callable.
1564 */
1565 template <typename IndexType>
1566 static void detect_super_cols
1567 (
1568   /* === Parameters ======================================================= */
1569 
1570   ColStructure<IndexType> Col [],    /* of size n_col+1 */
1571   IndexType A [],     /* row indices of A */
1572   IndexType head [],    /* head of degree lists and hash buckets */
1573   IndexType row_start,    /* pointer to set of columns to check */
1574   IndexType row_length    /* number of columns to check */
1575 )
1576 {
1577   /* === Local variables ================================================== */
1578 
1579   IndexType hash ;      /* hash value for a column */
1580   IndexType *rp ;     /* pointer to a row */
1581   IndexType c ;     /* a column index */
1582   IndexType super_c ;   /* column index of the column to absorb into */
1583   IndexType *cp1 ;      /* column pointer for column super_c */
1584   IndexType *cp2 ;      /* column pointer for column c */
1585   IndexType length ;    /* length of column super_c */
1586   IndexType prev_c ;    /* column preceding c in hash bucket */
1587   IndexType i ;     /* loop counter */
1588   IndexType *rp_end ;   /* pointer to the end of the row */
1589   IndexType col ;     /* a column index in the row to check */
1590   IndexType head_column ;   /* first column in hash bucket or degree list */
1591   IndexType first_col ;   /* first column in hash bucket */
1592 
1593   /* === Consider each column in the row ================================== */
1594 
1595   rp = &A [row_start] ;
1596   rp_end = rp + row_length ;
1597   while (rp < rp_end)
1598   {
1599     col = *rp++ ;
1600     if (Col[col].is_dead())
1601     {
1602       continue ;
1603     }
1604 
1605     /* get hash number for this column */
1606     hash = Col [col].shared3.hash ;
1607     COLAMD_ASSERT (hash <= n_col) ;
1608 
1609     /* === Get the first column in this hash bucket ===================== */
1610 
1611     head_column = head [hash] ;
1612     if (head_column > Empty)
1613     {
1614       first_col = Col [head_column].shared3.headhash ;
1615     }
1616     else
1617     {
1618       first_col = - (head_column + 2) ;
1619     }
1620 
1621     /* === Consider each column in the hash bucket ====================== */
1622 
1623     for (super_c = first_col ; super_c != Empty ;
1624      super_c = Col [super_c].shared4.hash_next)
1625     {
1626       COLAMD_ASSERT (Col [super_c].is_alive()) ;
1627       COLAMD_ASSERT (Col [super_c].shared3.hash == hash) ;
1628       length = Col [super_c].length ;
1629 
1630       /* prev_c is the column preceding column c in the hash bucket */
1631       prev_c = super_c ;
1632 
1633       /* === Compare super_c with all columns after it ================ */
1634 
1635       for (c = Col [super_c].shared4.hash_next ;
1636        c != Empty ; c = Col [c].shared4.hash_next)
1637       {
1638     COLAMD_ASSERT (c != super_c) ;
1639     COLAMD_ASSERT (Col[c].is_alive()) ;
1640     COLAMD_ASSERT (Col [c].shared3.hash == hash) ;
1641 
1642     /* not identical if lengths or scores are different */
1643     if (Col [c].length != length ||
1644         Col [c].shared2.score != Col [super_c].shared2.score)
1645     {
1646       prev_c = c ;
1647       continue ;
1648     }
1649 
1650     /* compare the two columns */
1651     cp1 = &A [Col [super_c].start] ;
1652     cp2 = &A [Col [c].start] ;
1653 
1654     for (i = 0 ; i < length ; i++)
1655     {
1656       /* the columns are "clean" (no dead rows) */
1657       COLAMD_ASSERT ( cp1->is_alive() );
1658       COLAMD_ASSERT ( cp2->is_alive() );
1659       /* row indices will same order for both supercols, */
1660       /* no gather scatter necessary */
1661       if (*cp1++ != *cp2++)
1662       {
1663         break ;
1664       }
1665     }
1666 
1667     /* the two columns are different if the for-loop "broke" */
1668     if (i != length)
1669     {
1670       prev_c = c ;
1671       continue ;
1672     }
1673 
1674     /* === Got it!  two columns are identical =================== */
1675 
1676     COLAMD_ASSERT (Col [c].shared2.score == Col [super_c].shared2.score) ;
1677 
1678     Col [super_c].shared1.thickness += Col [c].shared1.thickness ;
1679     Col [c].shared1.parent = super_c ;
1680     Col[c].kill_non_principal() ;
1681     /* order c later, in order_children() */
1682     Col [c].shared2.order = Empty ;
1683     /* remove c from hash bucket */
1684     Col [prev_c].shared4.hash_next = Col [c].shared4.hash_next ;
1685       }
1686     }
1687 
1688     /* === Empty this hash bucket ======================================= */
1689 
1690     if (head_column > Empty)
1691     {
1692       /* corresponding degree list "hash" is not empty */
1693       Col [head_column].shared3.headhash = Empty ;
1694     }
1695     else
1696     {
1697       /* corresponding degree list "hash" is empty */
1698       head [hash] = Empty ;
1699     }
1700   }
1701 }
1702 
1703 
1704 /* ========================================================================== */
1705 /* === garbage_collection =================================================== */
1706 /* ========================================================================== */
1707 
1708 /*
1709   Defragments and compacts columns and rows in the workspace A.  Used when
1710   all available memory has been used while performing row merging.  Returns
1711   the index of the first free position in A, after garbage collection.  The
1712   time taken by this routine is linear is the size of the array A, which is
1713   itself linear in the number of nonzeros in the input matrix.
1714   Not user-callable.
1715 */
1716 template <typename IndexType>
1717 static IndexType garbage_collection  /* returns the new value of pfree */
1718   (
1719     /* === Parameters ======================================================= */
1720 
1721     IndexType n_row,      /* number of rows */
1722     IndexType n_col,      /* number of columns */
1723     RowStructure<IndexType> Row [],    /* row info */
1724     ColStructure<IndexType> Col [],    /* column info */
1725     IndexType A [],     /* A [0 ... Alen-1] holds the matrix */
1726     IndexType *pfree      /* &A [0] ... pfree is in use */
1727     )
1728 {
1729   /* === Local variables ================================================== */
1730 
1731   IndexType *psrc ;     /* source pointer */
1732   IndexType *pdest ;    /* destination pointer */
1733   IndexType j ;     /* counter */
1734   IndexType r ;     /* a row index */
1735   IndexType c ;     /* a column index */
1736   IndexType length ;    /* length of a row or column */
1737 
1738   /* === Defragment the columns =========================================== */
1739 
1740   pdest = &A[0] ;
1741   for (c = 0 ; c < n_col ; c++)
1742   {
1743     if (Col[c].is_alive())
1744     {
1745       psrc = &A [Col [c].start] ;
1746 
1747       /* move and compact the column */
1748       COLAMD_ASSERT (pdest <= psrc) ;
1749       Col [c].start = (IndexType) (pdest - &A [0]) ;
1750       length = Col [c].length ;
1751       for (j = 0 ; j < length ; j++)
1752       {
1753     r = *psrc++ ;
1754     if (Row[r].is_alive())
1755     {
1756       *pdest++ = r ;
1757     }
1758       }
1759       Col [c].length = (IndexType) (pdest - &A [Col [c].start]) ;
1760     }
1761   }
1762 
1763   /* === Prepare to defragment the rows =================================== */
1764 
1765   for (r = 0 ; r < n_row ; r++)
1766   {
1767     if (Row[r].is_alive())
1768     {
1769       if (Row [r].length == 0)
1770       {
1771         /* this row is of zero length.  cannot compact it, so kill it */
1772         COLAMD_DEBUG3 (("Defrag row kill\n")) ;
1773         Row[r].kill() ;
1774       }
1775       else
1776       {
1777         /* save first column index in Row [r].shared2.first_column */
1778         psrc = &A [Row [r].start] ;
1779         Row [r].shared2.first_column = *psrc ;
1780         COLAMD_ASSERT (Row[r].is_alive()) ;
1781         /* flag the start of the row with the one's complement of row */
1782         *psrc = ones_complement(r) ;
1783 
1784       }
1785     }
1786   }
1787 
1788   /* === Defragment the rows ============================================== */
1789 
1790   psrc = pdest ;
1791   while (psrc < pfree)
1792   {
1793     /* find a negative number ... the start of a row */
1794     if (*psrc++ < 0)
1795     {
1796       psrc-- ;
1797       /* get the row index */
1798       r = ones_complement(*psrc) ;
1799       COLAMD_ASSERT (r >= 0 && r < n_row) ;
1800       /* restore first column index */
1801       *psrc = Row [r].shared2.first_column ;
1802       COLAMD_ASSERT (Row[r].is_alive()) ;
1803 
1804       /* move and compact the row */
1805       COLAMD_ASSERT (pdest <= psrc) ;
1806       Row [r].start = (IndexType) (pdest - &A [0]) ;
1807       length = Row [r].length ;
1808       for (j = 0 ; j < length ; j++)
1809       {
1810     c = *psrc++ ;
1811     if (Col[c].is_alive())
1812     {
1813       *pdest++ = c ;
1814     }
1815       }
1816       Row [r].length = (IndexType) (pdest - &A [Row [r].start]) ;
1817 
1818     }
1819   }
1820   /* ensure we found all the rows */
1821   COLAMD_ASSERT (debug_rows == 0) ;
1822 
1823   /* === Return the new value of pfree ==================================== */
1824 
1825   return ((IndexType) (pdest - &A [0])) ;
1826 }
1827 
1828 
1829 /* ========================================================================== */
1830 /* === clear_mark =========================================================== */
1831 /* ========================================================================== */
1832 
1833 /*
1834   Clears the Row [].shared2.mark array, and returns the new tag_mark.
1835   Return value is the new tag_mark.  Not user-callable.
1836 */
1837 template <typename IndexType>
1838 static inline  IndexType clear_mark  /* return the new value for tag_mark */
1839   (
1840       /* === Parameters ======================================================= */
1841 
1842     IndexType n_row,    /* number of rows in A */
1843     RowStructure<IndexType> Row [] /* Row [0 ... n_row-1].shared2.mark is set to zero */
1844     )
1845 {
1846   /* === Local variables ================================================== */
1847 
1848   IndexType r ;
1849 
1850   for (r = 0 ; r < n_row ; r++)
1851   {
1852     if (Row[r].is_alive())
1853     {
1854       Row [r].shared2.mark = 0 ;
1855     }
1856   }
1857   return (1) ;
1858 }
1859 
1860 } // namespace Colamd
1861 
1862 } // namespace internal
1863 #endif