Amd.h
Go to the documentation of this file.
00001 // This file is part of Eigen, a lightweight C++ template library
00002 // for linear algebra.
00003 //
00004 // Copyright (C) 2010 Gael Guennebaud <gael.guennebaud@inria.fr>
00005 
00006 /*
00007 
00008 NOTE: this routine has been adapted from the CSparse library:
00009 
00010 Copyright (c) 2006, Timothy A. Davis.
00011 http://www.cise.ufl.edu/research/sparse/CSparse
00012 
00013 CSparse is free software; you can redistribute it and/or
00014 modify it under the terms of the GNU Lesser General Public
00015 License as published by the Free Software Foundation; either
00016 version 2.1 of the License, or (at your option) any later version.
00017 
00018 CSparse is distributed in the hope that it will be useful,
00019 but WITHOUT ANY WARRANTY; without even the implied warranty of
00020 MERCHANTABILITY or FITNESS FOR A PARTICULAR PURPOSE.  See the GNU
00021 Lesser General Public License for more details.
00022 
00023 You should have received a copy of the GNU Lesser General Public
00024 License along with this Module; if not, write to the Free Software
00025 Foundation, Inc., 51 Franklin St, Fifth Floor, Boston, MA  02110-1301  USA
00026 
00027 */
00028 
00029 #include "../Core/util/NonMPL2.h"
00030 
00031 #ifndef EIGEN_SPARSE_AMD_H
00032 #define EIGEN_SPARSE_AMD_H
00033 
00034 namespace Eigen { 
00035 
00036 namespace internal {
00037   
00038 template<typename T> inline T amd_flip(const T& i) { return -i-2; }
00039 template<typename T> inline T amd_unflip(const T& i) { return i<0 ? amd_flip(i) : i; }
00040 template<typename T0, typename T1> inline bool amd_marked(const T0* w, const T1& j) { return w[j]<0; }
00041 template<typename T0, typename T1> inline void amd_mark(const T0* w, const T1& j) { return w[j] = amd_flip(w[j]); }
00042 
00043 /* clear w */
00044 template<typename Index>
00045 static int cs_wclear (Index mark, Index lemax, Index *w, Index n)
00046 {
00047   Index k;
00048   if(mark < 2 || (mark + lemax < 0))
00049   {
00050     for(k = 0; k < n; k++)
00051       if(w[k] != 0)
00052         w[k] = 1;
00053     mark = 2;
00054   }
00055   return (mark);     /* at this point, w[0..n-1] < mark holds */
00056 }
00057 
00058 /* depth-first search and postorder of a tree rooted at node j */
00059 template<typename Index>
00060 Index cs_tdfs(Index j, Index k, Index *head, const Index *next, Index *post, Index *stack)
00061 {
00062   int i, p, top = 0;
00063   if(!head || !next || !post || !stack) return (-1);    /* check inputs */
00064   stack[0] = j;                 /* place j on the stack */
00065   while (top >= 0)                /* while (stack is not empty) */
00066   {
00067     p = stack[top];           /* p = top of stack */
00068     i = head[p];              /* i = youngest child of p */
00069     if(i == -1)
00070     {
00071       top--;                 /* p has no unordered children left */
00072       post[k++] = p;        /* node p is the kth postordered node */
00073     }
00074     else
00075     {
00076       head[p] = next[i];   /* remove i from children of p */
00077       stack[++top] = i;     /* start dfs on child node i */
00078     }
00079   }
00080   return k;
00081 }
00082 
00083 
00090 template<typename Scalar, typename Index>
00091 void minimum_degree_ordering(SparseMatrix<Scalar,ColMajor,Index>& C, PermutationMatrix<Dynamic,Dynamic,Index>& perm)
00092 {
00093   using std::sqrt;
00094   
00095   int d, dk, dext, lemax = 0, e, elenk, eln, i, j, k, k1,
00096       k2, k3, jlast, ln, dense, nzmax, mindeg = 0, nvi, nvj, nvk, mark, wnvi,
00097       ok, nel = 0, p, p1, p2, p3, p4, pj, pk, pk1, pk2, pn, q, t;
00098   unsigned int h;
00099   
00100   Index n = C.cols();
00101   dense = std::max<Index> (16, Index(10 * sqrt(double(n))));   /* find dense threshold */
00102   dense = std::min<Index> (n-2, dense);
00103   
00104   Index cnz = C.nonZeros();
00105   perm.resize(n+1);
00106   t = cnz + cnz/5 + 2*n;                 /* add elbow room to C */
00107   C.resizeNonZeros(t);
00108   
00109   Index* W       = new Index[8*(n+1)]; /* get workspace */
00110   Index* len     = W;
00111   Index* nv      = W +   (n+1);
00112   Index* next    = W + 2*(n+1);
00113   Index* head    = W + 3*(n+1);
00114   Index* elen    = W + 4*(n+1);
00115   Index* degree  = W + 5*(n+1);
00116   Index* w       = W + 6*(n+1);
00117   Index* hhead   = W + 7*(n+1);
00118   Index* last    = perm.indices().data();                              /* use P as workspace for last */
00119   
00120   /* --- Initialize quotient graph ---------------------------------------- */
00121   Index* Cp = C.outerIndexPtr();
00122   Index* Ci = C.innerIndexPtr();
00123   for(k = 0; k < n; k++)
00124     len[k] = Cp[k+1] - Cp[k];
00125   len[n] = 0;
00126   nzmax = t;
00127   
00128   for(i = 0; i <= n; i++)
00129   {
00130     head[i]   = -1;                     // degree list i is empty
00131     last[i]   = -1;
00132     next[i]   = -1;
00133     hhead[i]  = -1;                     // hash list i is empty 
00134     nv[i]     = 1;                      // node i is just one node
00135     w[i]      = 1;                      // node i is alive
00136     elen[i]   = 0;                      // Ek of node i is empty
00137     degree[i] = len[i];                 // degree of node i
00138   }
00139   mark = internal::cs_wclear<Index>(0, 0, w, n);         /* clear w */
00140   elen[n] = -2;                         /* n is a dead element */
00141   Cp[n] = -1;                           /* n is a root of assembly tree */
00142   w[n] = 0;                             /* n is a dead element */
00143   
00144   /* --- Initialize degree lists ------------------------------------------ */
00145   for(i = 0; i < n; i++)
00146   {
00147     d = degree[i];
00148     if(d == 0)                         /* node i is empty */
00149     {
00150       elen[i] = -2;                 /* element i is dead */
00151       nel++;
00152       Cp[i] = -1;                   /* i is a root of assembly tree */
00153       w[i] = 0;
00154     }
00155     else if(d > dense)                 /* node i is dense */
00156     {
00157       nv[i] = 0;                    /* absorb i into element n */
00158       elen[i] = -1;                 /* node i is dead */
00159       nel++;
00160       Cp[i] = amd_flip (n);
00161       nv[n]++;
00162     }
00163     else
00164     {
00165       if(head[d] != -1) last[head[d]] = i;
00166       next[i] = head[d];           /* put node i in degree list d */
00167       head[d] = i;
00168     }
00169   }
00170   
00171   while (nel < n)                         /* while (selecting pivots) do */
00172   {
00173     /* --- Select node of minimum approximate degree -------------------- */
00174     for(k = -1; mindeg < n && (k = head[mindeg]) == -1; mindeg++) {}
00175     if(next[k] != -1) last[next[k]] = -1;
00176     head[mindeg] = next[k];          /* remove k from degree list */
00177     elenk = elen[k];                  /* elenk = |Ek| */
00178     nvk = nv[k];                      /* # of nodes k represents */
00179     nel += nvk;                        /* nv[k] nodes of A eliminated */
00180     
00181     /* --- Garbage collection ------------------------------------------- */
00182     if(elenk > 0 && cnz + mindeg >= nzmax)
00183     {
00184       for(j = 0; j < n; j++)
00185       {
00186         if((p = Cp[j]) >= 0)      /* j is a live node or element */
00187         {
00188           Cp[j] = Ci[p];          /* save first entry of object */
00189           Ci[p] = amd_flip (j);    /* first entry is now amd_flip(j) */
00190         }
00191       }
00192       for(q = 0, p = 0; p < cnz; ) /* scan all of memory */
00193       {
00194         if((j = amd_flip (Ci[p++])) >= 0)  /* found object j */
00195         {
00196           Ci[q] = Cp[j];       /* restore first entry of object */
00197           Cp[j] = q++;          /* new pointer to object j */
00198           for(k3 = 0; k3 < len[j]-1; k3++) Ci[q++] = Ci[p++];
00199         }
00200       }
00201       cnz = q;                       /* Ci[cnz...nzmax-1] now free */
00202     }
00203     
00204     /* --- Construct new element ---------------------------------------- */
00205     dk = 0;
00206     nv[k] = -nvk;                     /* flag k as in Lk */
00207     p = Cp[k];
00208     pk1 = (elenk == 0) ? p : cnz;      /* do in place if elen[k] == 0 */
00209     pk2 = pk1;
00210     for(k1 = 1; k1 <= elenk + 1; k1++)
00211     {
00212       if(k1 > elenk)
00213       {
00214         e = k;                     /* search the nodes in k */
00215         pj = p;                    /* list of nodes starts at Ci[pj]*/
00216         ln = len[k] - elenk;      /* length of list of nodes in k */
00217       }
00218       else
00219       {
00220         e = Ci[p++];              /* search the nodes in e */
00221         pj = Cp[e];
00222         ln = len[e];              /* length of list of nodes in e */
00223       }
00224       for(k2 = 1; k2 <= ln; k2++)
00225       {
00226         i = Ci[pj++];
00227         if((nvi = nv[i]) <= 0) continue; /* node i dead, or seen */
00228         dk += nvi;                 /* degree[Lk] += size of node i */
00229         nv[i] = -nvi;             /* negate nv[i] to denote i in Lk*/
00230         Ci[pk2++] = i;            /* place i in Lk */
00231         if(next[i] != -1) last[next[i]] = last[i];
00232         if(last[i] != -1)         /* remove i from degree list */
00233         {
00234           next[last[i]] = next[i];
00235         }
00236         else
00237         {
00238           head[degree[i]] = next[i];
00239         }
00240       }
00241       if(e != k)
00242       {
00243         Cp[e] = amd_flip (k);      /* absorb e into k */
00244         w[e] = 0;                 /* e is now a dead element */
00245       }
00246     }
00247     if(elenk != 0) cnz = pk2;         /* Ci[cnz...nzmax] is free */
00248     degree[k] = dk;                   /* external degree of k - |Lk\i| */
00249     Cp[k] = pk1;                      /* element k is in Ci[pk1..pk2-1] */
00250     len[k] = pk2 - pk1;
00251     elen[k] = -2;                     /* k is now an element */
00252     
00253     /* --- Find set differences ----------------------------------------- */
00254     mark = internal::cs_wclear<Index>(mark, lemax, w, n);  /* clear w if necessary */
00255     for(pk = pk1; pk < pk2; pk++)    /* scan 1: find |Le\Lk| */
00256     {
00257       i = Ci[pk];
00258       if((eln = elen[i]) <= 0) continue;/* skip if elen[i] empty */
00259       nvi = -nv[i];                      /* nv[i] was negated */
00260       wnvi = mark - nvi;
00261       for(p = Cp[i]; p <= Cp[i] + eln - 1; p++)  /* scan Ei */
00262       {
00263         e = Ci[p];
00264         if(w[e] >= mark)
00265         {
00266           w[e] -= nvi;          /* decrement |Le\Lk| */
00267         }
00268         else if(w[e] != 0)        /* ensure e is a live element */
00269         {
00270           w[e] = degree[e] + wnvi; /* 1st time e seen in scan 1 */
00271         }
00272       }
00273     }
00274     
00275     /* --- Degree update ------------------------------------------------ */
00276     for(pk = pk1; pk < pk2; pk++)    /* scan2: degree update */
00277     {
00278       i = Ci[pk];                   /* consider node i in Lk */
00279       p1 = Cp[i];
00280       p2 = p1 + elen[i] - 1;
00281       pn = p1;
00282       for(h = 0, d = 0, p = p1; p <= p2; p++)    /* scan Ei */
00283       {
00284         e = Ci[p];
00285         if(w[e] != 0)             /* e is an unabsorbed element */
00286         {
00287           dext = w[e] - mark;   /* dext = |Le\Lk| */
00288           if(dext > 0)
00289           {
00290             d += dext;         /* sum up the set differences */
00291             Ci[pn++] = e;     /* keep e in Ei */
00292             h += e;            /* compute the hash of node i */
00293           }
00294           else
00295           {
00296             Cp[e] = amd_flip (k);  /* aggressive absorb. e->k */
00297             w[e] = 0;             /* e is a dead element */
00298           }
00299         }
00300       }
00301       elen[i] = pn - p1 + 1;        /* elen[i] = |Ei| */
00302       p3 = pn;
00303       p4 = p1 + len[i];
00304       for(p = p2 + 1; p < p4; p++) /* prune edges in Ai */
00305       {
00306         j = Ci[p];
00307         if((nvj = nv[j]) <= 0) continue; /* node j dead or in Lk */
00308         d += nvj;                  /* degree(i) += |j| */
00309         Ci[pn++] = j;             /* place j in node list of i */
00310         h += j;                    /* compute hash for node i */
00311       }
00312       if(d == 0)                     /* check for mass elimination */
00313       {
00314         Cp[i] = amd_flip (k);      /* absorb i into k */
00315         nvi = -nv[i];
00316         dk -= nvi;                 /* |Lk| -= |i| */
00317         nvk += nvi;                /* |k| += nv[i] */
00318         nel += nvi;
00319         nv[i] = 0;
00320         elen[i] = -1;             /* node i is dead */
00321       }
00322       else
00323       {
00324         degree[i] = std::min<Index> (degree[i], d);   /* update degree(i) */
00325         Ci[pn] = Ci[p3];         /* move first node to end */
00326         Ci[p3] = Ci[p1];         /* move 1st el. to end of Ei */
00327         Ci[p1] = k;               /* add k as 1st element in of Ei */
00328         len[i] = pn - p1 + 1;     /* new len of adj. list of node i */
00329         h %= n;                    /* finalize hash of i */
00330         next[i] = hhead[h];      /* place i in hash bucket */
00331         hhead[h] = i;
00332         last[i] = h;              /* save hash of i in last[i] */
00333       }
00334     }                                   /* scan2 is done */
00335     degree[k] = dk;                   /* finalize |Lk| */
00336     lemax = std::max<Index>(lemax, dk);
00337     mark = internal::cs_wclear<Index>(mark+lemax, lemax, w, n);    /* clear w */
00338     
00339     /* --- Supernode detection ------------------------------------------ */
00340     for(pk = pk1; pk < pk2; pk++)
00341     {
00342       i = Ci[pk];
00343       if(nv[i] >= 0) continue;         /* skip if i is dead */
00344       h = last[i];                      /* scan hash bucket of node i */
00345       i = hhead[h];
00346       hhead[h] = -1;                    /* hash bucket will be empty */
00347       for(; i != -1 && next[i] != -1; i = next[i], mark++)
00348       {
00349         ln = len[i];
00350         eln = elen[i];
00351         for(p = Cp[i]+1; p <= Cp[i] + ln-1; p++) w[Ci[p]] = mark;
00352         jlast = i;
00353         for(j = next[i]; j != -1; ) /* compare i with all j */
00354         {
00355           ok = (len[j] == ln) && (elen[j] == eln);
00356           for(p = Cp[j] + 1; ok && p <= Cp[j] + ln - 1; p++)
00357           {
00358             if(w[Ci[p]] != mark) ok = 0;    /* compare i and j*/
00359           }
00360           if(ok)                     /* i and j are identical */
00361           {
00362             Cp[j] = amd_flip (i);  /* absorb j into i */
00363             nv[i] += nv[j];
00364             nv[j] = 0;
00365             elen[j] = -1;         /* node j is dead */
00366             j = next[j];          /* delete j from hash bucket */
00367             next[jlast] = j;
00368           }
00369           else
00370           {
00371             jlast = j;             /* j and i are different */
00372             j = next[j];
00373           }
00374         }
00375       }
00376     }
00377     
00378     /* --- Finalize new element------------------------------------------ */
00379     for(p = pk1, pk = pk1; pk < pk2; pk++)   /* finalize Lk */
00380     {
00381       i = Ci[pk];
00382       if((nvi = -nv[i]) <= 0) continue;/* skip if i is dead */
00383       nv[i] = nvi;                      /* restore nv[i] */
00384       d = degree[i] + dk - nvi;         /* compute external degree(i) */
00385       d = std::min<Index> (d, n - nel - nvi);
00386       if(head[d] != -1) last[head[d]] = i;
00387       next[i] = head[d];               /* put i back in degree list */
00388       last[i] = -1;
00389       head[d] = i;
00390       mindeg = std::min<Index> (mindeg, d);       /* find new minimum degree */
00391       degree[i] = d;
00392       Ci[p++] = i;                      /* place i in Lk */
00393     }
00394     nv[k] = nvk;                      /* # nodes absorbed into k */
00395     if((len[k] = p-pk1) == 0)         /* length of adj list of element k*/
00396     {
00397       Cp[k] = -1;                   /* k is a root of the tree */
00398       w[k] = 0;                     /* k is now a dead element */
00399     }
00400     if(elenk != 0) cnz = p;           /* free unused space in Lk */
00401   }
00402   
00403   /* --- Postordering ----------------------------------------------------- */
00404   for(i = 0; i < n; i++) Cp[i] = amd_flip (Cp[i]);/* fix assembly tree */
00405   for(j = 0; j <= n; j++) head[j] = -1;
00406   for(j = n; j >= 0; j--)              /* place unordered nodes in lists */
00407   {
00408     if(nv[j] > 0) continue;          /* skip if j is an element */
00409     next[j] = head[Cp[j]];          /* place j in list of its parent */
00410     head[Cp[j]] = j;
00411   }
00412   for(e = n; e >= 0; e--)              /* place elements in lists */
00413   {
00414     if(nv[e] <= 0) continue;         /* skip unless e is an element */
00415     if(Cp[e] != -1)
00416     {
00417       next[e] = head[Cp[e]];      /* place e in list of its parent */
00418       head[Cp[e]] = e;
00419     }
00420   }
00421   for(k = 0, i = 0; i <= n; i++)       /* postorder the assembly tree */
00422   {
00423     if(Cp[i] == -1) k = internal::cs_tdfs<Index>(i, k, head, next, perm.indices().data(), w);
00424   }
00425   
00426   perm.indices().conservativeResize(n);
00427 
00428   delete[] W;
00429 }
00430 
00431 } // namespace internal
00432 
00433 } // end namespace Eigen
00434 
00435 #endif // EIGEN_SPARSE_AMD_H


acado
Author(s): Milan Vukov, Rien Quirynen
autogenerated on Thu Aug 27 2015 11:57:48