Eigen  3.4.90 (git rev 5a9f66fb35d03a4da9ef8976e67a61b30aa16dcf)
 
Loading...
Searching...
No Matches
Amd.h
1// This file is part of Eigen, a lightweight C++ template library
2// for linear algebra.
3//
4// Copyright (C) 2010 Gael Guennebaud <[email protected]>
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/*
11NOTE: this routine has been adapted from the CSparse library:
12
13Copyright (c) 2006, Timothy A. Davis.
14http://www.suitesparse.com
15
16The author of CSparse, Timothy A. Davis., has executed a license with Google LLC
17to permit distribution of this code and derivative works as part of Eigen under
18the Mozilla Public License v. 2.0, as stated at the top of this file.
19*/
20
21#ifndef EIGEN_SPARSE_AMD_H
22#define EIGEN_SPARSE_AMD_H
23
24// IWYU pragma: private
25#include "./InternalHeaderCheck.h"
26
27namespace Eigen {
28
29namespace internal {
30
31template <typename T>
32inline T amd_flip(const T& i) {
33 return -i - 2;
34}
35template <typename T>
36inline T amd_unflip(const T& i) {
37 return i < 0 ? amd_flip(i) : i;
38}
39template <typename T0, typename T1>
40inline bool amd_marked(const T0* w, const T1& j) {
41 return w[j] < 0;
42}
43template <typename T0, typename T1>
44inline void amd_mark(const T0* w, const T1& j) {
45 return w[j] = amd_flip(w[j]);
46}
47
48/* clear w */
49template <typename StorageIndex>
50static StorageIndex cs_wclear(StorageIndex mark, StorageIndex lemax, StorageIndex* w, StorageIndex n) {
51 StorageIndex k;
52 if (mark < 2 || (mark + lemax < 0)) {
53 for (k = 0; k < n; k++)
54 if (w[k] != 0) w[k] = 1;
55 mark = 2;
56 }
57 return (mark); /* at this point, w[0..n-1] < mark holds */
58}
59
60/* depth-first search and postorder of a tree rooted at node j */
61template <typename StorageIndex>
62StorageIndex cs_tdfs(StorageIndex j, StorageIndex k, StorageIndex* head, const StorageIndex* next, StorageIndex* post,
63 StorageIndex* stack) {
64 StorageIndex i, p, top = 0;
65 if (!head || !next || !post || !stack) return (-1); /* check inputs */
66 stack[0] = j; /* place j on the stack */
67 while (top >= 0) /* while (stack is not empty) */
68 {
69 p = stack[top]; /* p = top of stack */
70 i = head[p]; /* i = youngest child of p */
71 if (i == -1) {
72 top--; /* p has no unordered children left */
73 post[k++] = p; /* node p is the kth postordered node */
74 } else {
75 head[p] = next[i]; /* remove i from children of p */
76 stack[++top] = i; /* start dfs on child node i */
77 }
78 }
79 return k;
80}
81
91template <typename Scalar, typename StorageIndex>
92void minimum_degree_ordering(SparseMatrix<Scalar, ColMajor, StorageIndex>& C,
93 PermutationMatrix<Dynamic, Dynamic, StorageIndex>& perm) {
94 using std::sqrt;
95
96 StorageIndex d, dk, dext, lemax = 0, e, elenk, eln, i, j, k, k1, k2, k3, jlast, ln, dense, nzmax, mindeg = 0, nvi,
97 nvj, nvk, mark, wnvi, ok, nel = 0, p, p1, p2, p3, p4, pj, pk, pk1, pk2, pn, q, t, h;
98
99 StorageIndex n = StorageIndex(C.cols());
100 dense = std::max<StorageIndex>(16, StorageIndex(10 * sqrt(double(n)))); /* find dense threshold */
101 dense = (std::min)(n - 2, dense);
102
103 StorageIndex cnz = StorageIndex(C.nonZeros());
104 perm.resize(n + 1);
105 t = cnz + cnz / 5 + 2 * n; /* add elbow room to C */
106 C.resizeNonZeros(t);
107
108 // get workspace
109 ei_declare_aligned_stack_constructed_variable(StorageIndex, W, 8 * (n + 1), 0);
110 StorageIndex* len = W;
111 StorageIndex* nv = W + (n + 1);
112 StorageIndex* next = W + 2 * (n + 1);
113 StorageIndex* head = W + 3 * (n + 1);
114 StorageIndex* elen = W + 4 * (n + 1);
115 StorageIndex* degree = W + 5 * (n + 1);
116 StorageIndex* w = W + 6 * (n + 1);
117 StorageIndex* hhead = W + 7 * (n + 1);
118 StorageIndex* last = perm.indices().data(); /* use P as workspace for last */
119
120 /* --- Initialize quotient graph ---------------------------------------- */
121 StorageIndex* Cp = C.outerIndexPtr();
122 StorageIndex* Ci = C.innerIndexPtr();
123 for (k = 0; k < n; k++) len[k] = Cp[k + 1] - Cp[k];
124 len[n] = 0;
125 nzmax = t;
126
127 for (i = 0; i <= n; i++) {
128 head[i] = -1; // degree list i is empty
129 last[i] = -1;
130 next[i] = -1;
131 hhead[i] = -1; // hash list i is empty
132 nv[i] = 1; // node i is just one node
133 w[i] = 1; // node i is alive
134 elen[i] = 0; // Ek of node i is empty
135 degree[i] = len[i]; // degree of node i
136 }
137 mark = internal::cs_wclear<StorageIndex>(0, 0, w, n); /* clear w */
138
139 /* --- Initialize degree lists ------------------------------------------ */
140 for (i = 0; i < n; i++) {
141 bool has_diag = false;
142 for (p = Cp[i]; p < Cp[i + 1]; ++p)
143 if (Ci[p] == i) {
144 has_diag = true;
145 break;
146 }
147
148 d = degree[i];
149 if (d == 1 && has_diag) /* node i is empty */
150 {
151 elen[i] = -2; /* element i is dead */
152 nel++;
153 Cp[i] = -1; /* i is a root of assembly tree */
154 w[i] = 0;
155 } else if (d > dense || !has_diag) /* node i is dense or has no structural diagonal element */
156 {
157 nv[i] = 0; /* absorb i into element n */
158 elen[i] = -1; /* node i is dead */
159 nel++;
160 Cp[i] = amd_flip(n);
161 nv[n]++;
162 } else {
163 if (head[d] != -1) last[head[d]] = i;
164 next[i] = head[d]; /* put node i in degree list d */
165 head[d] = i;
166 }
167 }
168
169 elen[n] = -2; /* n is a dead element */
170 Cp[n] = -1; /* n is a root of assembly tree */
171 w[n] = 0; /* n is a dead element */
172
173 while (nel < n) /* while (selecting pivots) do */
174 {
175 /* --- Select node of minimum approximate degree -------------------- */
176 for (k = -1; mindeg < n && (k = head[mindeg]) == -1; mindeg++) {
177 }
178 if (next[k] != -1) last[next[k]] = -1;
179 head[mindeg] = next[k]; /* remove k from degree list */
180 elenk = elen[k]; /* elenk = |Ek| */
181 nvk = nv[k]; /* # of nodes k represents */
182 nel += nvk; /* nv[k] nodes of A eliminated */
183
184 /* --- Garbage collection ------------------------------------------- */
185 if (elenk > 0 && cnz + mindeg >= nzmax) {
186 for (j = 0; j < n; j++) {
187 if ((p = Cp[j]) >= 0) /* j is a live node or element */
188 {
189 Cp[j] = Ci[p]; /* save first entry of object */
190 Ci[p] = amd_flip(j); /* first entry is now amd_flip(j) */
191 }
192 }
193 for (q = 0, p = 0; p < cnz;) /* scan all of memory */
194 {
195 if ((j = amd_flip(Ci[p++])) >= 0) /* found object j */
196 {
197 Ci[q] = Cp[j]; /* restore first entry of object */
198 Cp[j] = q++; /* new pointer to object j */
199 for (k3 = 0; k3 < len[j] - 1; k3++) Ci[q++] = Ci[p++];
200 }
201 }
202 cnz = q; /* Ci[cnz...nzmax-1] now free */
203 }
204
205 /* --- Construct new element ---------------------------------------- */
206 dk = 0;
207 nv[k] = -nvk; /* flag k as in Lk */
208 p = Cp[k];
209 pk1 = (elenk == 0) ? p : cnz; /* do in place if elen[k] == 0 */
210 pk2 = pk1;
211 for (k1 = 1; k1 <= elenk + 1; k1++) {
212 if (k1 > elenk) {
213 e = k; /* search the nodes in k */
214 pj = p; /* list of nodes starts at Ci[pj]*/
215 ln = len[k] - elenk; /* length of list of nodes in k */
216 } else {
217 e = Ci[p++]; /* search the nodes in e */
218 pj = Cp[e];
219 ln = len[e]; /* length of list of nodes in e */
220 }
221 for (k2 = 1; k2 <= ln; k2++) {
222 i = Ci[pj++];
223 if ((nvi = nv[i]) <= 0) continue; /* node i dead, or seen */
224 dk += nvi; /* degree[Lk] += size of node i */
225 nv[i] = -nvi; /* negate nv[i] to denote i in Lk*/
226 Ci[pk2++] = i; /* place i in Lk */
227 if (next[i] != -1) last[next[i]] = last[i];
228 if (last[i] != -1) /* remove i from degree list */
229 {
230 next[last[i]] = next[i];
231 } else {
232 head[degree[i]] = next[i];
233 }
234 }
235 if (e != k) {
236 Cp[e] = amd_flip(k); /* absorb e into k */
237 w[e] = 0; /* e is now a dead element */
238 }
239 }
240 if (elenk != 0) cnz = pk2; /* Ci[cnz...nzmax] is free */
241 degree[k] = dk; /* external degree of k - |Lk\i| */
242 Cp[k] = pk1; /* element k is in Ci[pk1..pk2-1] */
243 len[k] = pk2 - pk1;
244 elen[k] = -2; /* k is now an element */
245
246 /* --- Find set differences ----------------------------------------- */
247 mark = internal::cs_wclear<StorageIndex>(mark, lemax, w, n); /* clear w if necessary */
248 for (pk = pk1; pk < pk2; pk++) /* scan 1: find |Le\Lk| */
249 {
250 i = Ci[pk];
251 if ((eln = elen[i]) <= 0) continue; /* skip if elen[i] empty */
252 nvi = -nv[i]; /* nv[i] was negated */
253 wnvi = mark - nvi;
254 for (p = Cp[i]; p <= Cp[i] + eln - 1; p++) /* scan Ei */
255 {
256 e = Ci[p];
257 if (w[e] >= mark) {
258 w[e] -= nvi; /* decrement |Le\Lk| */
259 } else if (w[e] != 0) /* ensure e is a live element */
260 {
261 w[e] = degree[e] + wnvi; /* 1st time e seen in scan 1 */
262 }
263 }
264 }
265
266 /* --- Degree update ------------------------------------------------ */
267 for (pk = pk1; pk < pk2; pk++) /* scan2: degree update */
268 {
269 i = Ci[pk]; /* consider node i in Lk */
270 p1 = Cp[i];
271 p2 = p1 + elen[i] - 1;
272 pn = p1;
273 for (h = 0, d = 0, p = p1; p <= p2; p++) /* scan Ei */
274 {
275 e = Ci[p];
276 if (w[e] != 0) /* e is an unabsorbed element */
277 {
278 dext = w[e] - mark; /* dext = |Le\Lk| */
279 if (dext > 0) {
280 d += dext; /* sum up the set differences */
281 Ci[pn++] = e; /* keep e in Ei */
282 h += e; /* compute the hash of node i */
283 } else {
284 Cp[e] = amd_flip(k); /* aggressive absorb. e->k */
285 w[e] = 0; /* e is a dead element */
286 }
287 }
288 }
289 elen[i] = pn - p1 + 1; /* elen[i] = |Ei| */
290 p3 = pn;
291 p4 = p1 + len[i];
292 for (p = p2 + 1; p < p4; p++) /* prune edges in Ai */
293 {
294 j = Ci[p];
295 if ((nvj = nv[j]) <= 0) continue; /* node j dead or in Lk */
296 d += nvj; /* degree(i) += |j| */
297 Ci[pn++] = j; /* place j in node list of i */
298 h += j; /* compute hash for node i */
299 }
300 if (d == 0) /* check for mass elimination */
301 {
302 Cp[i] = amd_flip(k); /* absorb i into k */
303 nvi = -nv[i];
304 dk -= nvi; /* |Lk| -= |i| */
305 nvk += nvi; /* |k| += nv[i] */
306 nel += nvi;
307 nv[i] = 0;
308 elen[i] = -1; /* node i is dead */
309 } else {
310 degree[i] = std::min<StorageIndex>(degree[i], d); /* update degree(i) */
311 Ci[pn] = Ci[p3]; /* move first node to end */
312 Ci[p3] = Ci[p1]; /* move 1st el. to end of Ei */
313 Ci[p1] = k; /* add k as 1st element in of Ei */
314 len[i] = pn - p1 + 1; /* new len of adj. list of node i */
315 h %= n; /* finalize hash of i */
316 next[i] = hhead[h]; /* place i in hash bucket */
317 hhead[h] = i;
318 last[i] = h; /* save hash of i in last[i] */
319 }
320 } /* scan2 is done */
321 degree[k] = dk; /* finalize |Lk| */
322 lemax = std::max<StorageIndex>(lemax, dk);
323 mark = internal::cs_wclear<StorageIndex>(mark + lemax, lemax, w, n); /* clear w */
324
325 /* --- Supernode detection ------------------------------------------ */
326 for (pk = pk1; pk < pk2; pk++) {
327 i = Ci[pk];
328 if (nv[i] >= 0) continue; /* skip if i is dead */
329 h = last[i]; /* scan hash bucket of node i */
330 i = hhead[h];
331 hhead[h] = -1; /* hash bucket will be empty */
332 for (; i != -1 && next[i] != -1; i = next[i], mark++) {
333 ln = len[i];
334 eln = elen[i];
335 for (p = Cp[i] + 1; p <= Cp[i] + ln - 1; p++) w[Ci[p]] = mark;
336 jlast = i;
337 for (j = next[i]; j != -1;) /* compare i with all j */
338 {
339 ok = (len[j] == ln) && (elen[j] == eln);
340 for (p = Cp[j] + 1; ok && p <= Cp[j] + ln - 1; p++) {
341 if (w[Ci[p]] != mark) ok = 0; /* compare i and j*/
342 }
343 if (ok) /* i and j are identical */
344 {
345 Cp[j] = amd_flip(i); /* absorb j into i */
346 nv[i] += nv[j];
347 nv[j] = 0;
348 elen[j] = -1; /* node j is dead */
349 j = next[j]; /* delete j from hash bucket */
350 next[jlast] = j;
351 } else {
352 jlast = j; /* j and i are different */
353 j = next[j];
354 }
355 }
356 }
357 }
358
359 /* --- Finalize new element------------------------------------------ */
360 for (p = pk1, pk = pk1; pk < pk2; pk++) /* finalize Lk */
361 {
362 i = Ci[pk];
363 if ((nvi = -nv[i]) <= 0) continue; /* skip if i is dead */
364 nv[i] = nvi; /* restore nv[i] */
365 d = degree[i] + dk - nvi; /* compute external degree(i) */
366 d = std::min<StorageIndex>(d, n - nel - nvi);
367 if (head[d] != -1) last[head[d]] = i;
368 next[i] = head[d]; /* put i back in degree list */
369 last[i] = -1;
370 head[d] = i;
371 mindeg = std::min<StorageIndex>(mindeg, d); /* find new minimum degree */
372 degree[i] = d;
373 Ci[p++] = i; /* place i in Lk */
374 }
375 nv[k] = nvk; /* # nodes absorbed into k */
376 if ((len[k] = p - pk1) == 0) /* length of adj list of element k*/
377 {
378 Cp[k] = -1; /* k is a root of the tree */
379 w[k] = 0; /* k is now a dead element */
380 }
381 if (elenk != 0) cnz = p; /* free unused space in Lk */
382 }
383
384 /* --- Postordering ----------------------------------------------------- */
385 for (i = 0; i < n; i++) Cp[i] = amd_flip(Cp[i]); /* fix assembly tree */
386 for (j = 0; j <= n; j++) head[j] = -1;
387 for (j = n; j >= 0; j--) /* place unordered nodes in lists */
388 {
389 if (nv[j] > 0) continue; /* skip if j is an element */
390 next[j] = head[Cp[j]]; /* place j in list of its parent */
391 head[Cp[j]] = j;
392 }
393 for (e = n; e >= 0; e--) /* place elements in lists */
394 {
395 if (nv[e] <= 0) continue; /* skip unless e is an element */
396 if (Cp[e] != -1) {
397 next[e] = head[Cp[e]]; /* place e in list of its parent */
398 head[Cp[e]] = e;
399 }
400 }
401 for (k = 0, i = 0; i <= n; i++) /* postorder the assembly tree */
402 {
403 if (Cp[i] == -1) k = internal::cs_tdfs<StorageIndex>(i, k, head, next, perm.indices().data(), w);
404 }
405
406 perm.indices().conservativeResize(n);
407}
408
409} // namespace internal
410
411} // end namespace Eigen
412
413#endif // EIGEN_SPARSE_AMD_H
Namespace containing all symbols from the Eigen library.
Definition Core:137
const Eigen::CwiseUnaryOp< Eigen::internal::scalar_sqrt_op< typename Derived::Scalar >, const Derived > sqrt(const Eigen::ArrayBase< Derived > &x)