Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
amesos_amd.h
Go to the documentation of this file.
1 /* ========================================================================= */
2 /* === AMD: approximate minimum degree ordering =========================== */
3 /* ========================================================================= */
4 
5 /* ------------------------------------------------------------------------- */
6 /* AMD Version 2.2, Copyright (c) 2007 by Timothy A. Davis, */
7 /* Patrick R. Amestoy, and Iain S. Duff. See ../README.txt for License. */
8 /* email: davis at cise.ufl.edu CISE Department, Univ. of Florida. */
9 /* web: http://www.cise.ufl.edu/research/sparse/amd */
10 /* ------------------------------------------------------------------------- */
11 
12 /* AMD finds a symmetric ordering P of a matrix A so that the Cholesky
13  * factorization of P*A*P' has fewer nonzeros and takes less work than the
14  * Cholesky factorization of A. If A is not symmetric, then it performs its
15  * ordering on the matrix A+A'. Two sets of user-callable routines are
16  * provided, one for int integers and the other for UF_long integers.
17  *
18  * The method is based on the approximate minimum degree algorithm, discussed
19  * in Amestoy, Davis, and Duff, "An approximate degree ordering algorithm",
20  * SIAM Journal of Matrix Analysis and Applications, vol. 17, no. 4, pp.
21  * 886-905, 1996. This package can perform both the AMD ordering (with
22  * aggressive absorption), and the AMDBAR ordering (without aggressive
23  * absorption) discussed in the above paper. This package differs from the
24  * Fortran codes discussed in the paper:
25  *
26  * (1) it can ignore "dense" rows and columns, leading to faster run times
27  * (2) it computes the ordering of A+A' if A is not symmetric
28  * (3) it is followed by a depth-first post-ordering of the assembly tree
29  * (or supernodal elimination tree)
30  *
31  * For historical reasons, the Fortran versions, amd.f and amdbar.f, have
32  * been left (nearly) unchanged. They compute the identical ordering as
33  * described in the above paper.
34  */
35 
36 #ifndef AMESOS_AMD_H
37 #define AMESOS_AMD_H
38 
39 /* make it easy for C++ programs to include AMD */
40 #ifdef __cplusplus
41 extern "C" {
42 #endif
43 
44 /* get the definition of size_t: */
45 #include <stddef.h>
46 
47 /* define UF_long */
48 #include "amesos_UFconfig.h"
49 
50 int amesos_amd_order /* returns AMD_OK, AMD_OK_BUT_JUMBLED,
51  * AMD_INVALID, or AMD_OUT_OF_MEMORY */
52 (
53  int n, /* A is n-by-n. n must be >= 0. */
54  const int Ap [ ], /* column pointers for A, of size n+1 */
55  const int Ai [ ], /* row indices of A, of size nz = Ap [n] */
56  int P [ ], /* output permutation, of size n */
57  double Control [ ], /* input Control settings, of size AMD_CONTROL */
58  double Info [ ] /* output Info statistics, of size AMD_INFO */
59 ) ;
60 
61 UF_long amesos_amd_l_order /* see above for description of arguments */
62 (
63  UF_long n,
64  const UF_long Ap [ ],
65  const UF_long Ai [ ],
66  UF_long P [ ],
67  double Control [ ],
68  double Info [ ]
69 ) ;
70 
71 /* Input arguments (not modified):
72  *
73  * n: the matrix A is n-by-n.
74  * Ap: an int/UF_long array of size n+1, containing column pointers of A.
75  * Ai: an int/UF_long array of size nz, containing the row indices of A,
76  * where nz = Ap [n].
77  * Control: a double array of size AMD_CONTROL, containing control
78  * parameters. Defaults are used if Control is NULL.
79  *
80  * Output arguments (not defined on input):
81  *
82  * P: an int/UF_long array of size n, containing the output permutation. If
83  * row i is the kth pivot row, then P [k] = i. In MATLAB notation,
84  * the reordered matrix is A (P,P).
85  * Info: a double array of size AMD_INFO, containing statistical
86  * information. Ignored if Info is NULL.
87  *
88  * On input, the matrix A is stored in column-oriented form. The row indices
89  * of nonzero entries in column j are stored in Ai [Ap [j] ... Ap [j+1]-1].
90  *
91  * If the row indices appear in ascending order in each column, and there
92  * are no duplicate entries, then amesos_amd_order is slightly more efficient in
93  * terms of time and memory usage. If this condition does not hold, a copy
94  * of the matrix is created (where these conditions do hold), and the copy is
95  * ordered. This feature is new to v2.0 (v1.2 and earlier required this
96  * condition to hold for the input matrix).
97  *
98  * Row indices must be in the range 0 to
99  * n-1. Ap [0] must be zero, and thus nz = Ap [n] is the number of nonzeros
100  * in A. The array Ap is of size n+1, and the array Ai is of size nz = Ap [n].
101  * The matrix does not need to be symmetric, and the diagonal does not need to
102  * be present (if diagonal entries are present, they are ignored except for
103  * the output statistic Info [AMD_NZDIAG]). The arrays Ai and Ap are not
104  * modified. This form of the Ap and Ai arrays to represent the nonzero
105  * pattern of the matrix A is the same as that used internally by MATLAB.
106  * If you wish to use a more flexible input structure, please see the
107  * umfpack_*_triplet_to_col routines in the UMFPACK package, at
108  * http://www.cise.ufl.edu/research/sparse/umfpack.
109  *
110  * Restrictions: n >= 0. Ap [0] = 0. Ap [j] <= Ap [j+1] for all j in the
111  * range 0 to n-1. nz = Ap [n] >= 0. Ai [0..nz-1] must be in the range 0
112  * to n-1. Finally, Ai, Ap, and P must not be NULL. If any of these
113  * restrictions are not met, AMD returns AMD_INVALID.
114  *
115  * AMD returns:
116  *
117  * AMD_OK if the matrix is valid and sufficient memory can be allocated to
118  * perform the ordering.
119  *
120  * AMD_OUT_OF_MEMORY if not enough memory can be allocated.
121  *
122  * AMD_INVALID if the input arguments n, Ap, Ai are invalid, or if P is
123  * NULL.
124  *
125  * AMD_OK_BUT_JUMBLED if the matrix had unsorted columns, and/or duplicate
126  * entries, but was otherwise valid.
127  *
128  * The AMD routine first forms the pattern of the matrix A+A', and then
129  * computes a fill-reducing ordering, P. If P [k] = i, then row/column i of
130  * the original is the kth pivotal row. In MATLAB notation, the permuted
131  * matrix is A (P,P), except that 0-based indexing is used instead of the
132  * 1-based indexing in MATLAB.
133  *
134  * The Control array is used to set various parameters for AMD. If a NULL
135  * pointer is passed, default values are used. The Control array is not
136  * modified.
137  *
138  * Control [AMD_DENSE]: controls the threshold for "dense" rows/columns.
139  * A dense row/column in A+A' can cause AMD to spend a lot of time in
140  * ordering the matrix. If Control [AMD_DENSE] >= 0, rows/columns
141  * with more than Control [AMD_DENSE] * sqrt (n) entries are ignored
142  * during the ordering, and placed last in the output order. The
143  * default value of Control [AMD_DENSE] is 10. If negative, no
144  * rows/columns are treated as "dense". Rows/columns with 16 or
145  * fewer off-diagonal entries are never considered "dense".
146  *
147  * Control [AMD_AGGRESSIVE]: controls whether or not to use aggressive
148  * absorption, in which a prior element is absorbed into the current
149  * element if is a subset of the current element, even if it is not
150  * adjacent to the current pivot element (refer to Amestoy, Davis,
151  * & Duff, 1996, for more details). The default value is nonzero,
152  * which means to perform aggressive absorption. This nearly always
153  * leads to a better ordering (because the approximate degrees are
154  * more accurate) and a lower execution time. There are cases where
155  * it can lead to a slightly worse ordering, however. To turn it off,
156  * set Control [AMD_AGGRESSIVE] to 0.
157  *
158  * Control [2..4] are not used in the current version, but may be used in
159  * future versions.
160  *
161  * The Info array provides statistics about the ordering on output. If it is
162  * not present, the statistics are not returned. This is not an error
163  * condition.
164  *
165  * Info [AMD_STATUS]: the return value of AMD, either AMD_OK,
166  * AMD_OK_BUT_JUMBLED, AMD_OUT_OF_MEMORY, or AMD_INVALID.
167  *
168  * Info [AMD_N]: n, the size of the input matrix
169  *
170  * Info [AMD_NZ]: the number of nonzeros in A, nz = Ap [n]
171  *
172  * Info [AMD_SYMMETRY]: the symmetry of the matrix A. It is the number
173  * of "matched" off-diagonal entries divided by the total number of
174  * off-diagonal entries. An entry A(i,j) is matched if A(j,i) is also
175  * an entry, for any pair (i,j) for which i != j. In MATLAB notation,
176  * S = spones (A) ;
177  * B = tril (S, -1) + triu (S, 1) ;
178  * symmetry = nnz (B & B') / nnz (B) ;
179  *
180  * Info [AMD_NZDIAG]: the number of entries on the diagonal of A.
181  *
182  * Info [AMD_NZ_A_PLUS_AT]: the number of nonzeros in A+A', excluding the
183  * diagonal. If A is perfectly symmetric (Info [AMD_SYMMETRY] = 1)
184  * with a fully nonzero diagonal, then Info [AMD_NZ_A_PLUS_AT] = nz-n
185  * (the smallest possible value). If A is perfectly unsymmetric
186  * (Info [AMD_SYMMETRY] = 0, for an upper triangular matrix, for
187  * example) with no diagonal, then Info [AMD_NZ_A_PLUS_AT] = 2*nz
188  * (the largest possible value).
189  *
190  * Info [AMD_NDENSE]: the number of "dense" rows/columns of A+A' that were
191  * removed from A prior to ordering. These are placed last in the
192  * output order P.
193  *
194  * Info [AMD_MEMORY]: the amount of memory used by AMD, in bytes. In the
195  * current version, this is 1.2 * Info [AMD_NZ_A_PLUS_AT] + 9*n
196  * times the size of an integer. This is at most 2.4nz + 9n. This
197  * excludes the size of the input arguments Ai, Ap, and P, which have
198  * a total size of nz + 2*n + 1 integers.
199  *
200  * Info [AMD_NCMPA]: the number of garbage collections performed.
201  *
202  * Info [AMD_LNZ]: the number of nonzeros in L (excluding the diagonal).
203  * This is a slight upper bound because mass elimination is combined
204  * with the approximate degree update. It is a rough upper bound if
205  * there are many "dense" rows/columns. The rest of the statistics,
206  * below, are also slight or rough upper bounds, for the same reasons.
207  * The post-ordering of the assembly tree might also not exactly
208  * correspond to a true elimination tree postordering.
209  *
210  * Info [AMD_NDIV]: the number of divide operations for a subsequent LDL'
211  * or LU factorization of the permuted matrix A (P,P).
212  *
213  * Info [AMD_NMULTSUBS_LDL]: the number of multiply-subtract pairs for a
214  * subsequent LDL' factorization of A (P,P).
215  *
216  * Info [AMD_NMULTSUBS_LU]: the number of multiply-subtract pairs for a
217  * subsequent LU factorization of A (P,P), assuming that no numerical
218  * pivoting is required.
219  *
220  * Info [AMD_DMAX]: the maximum number of nonzeros in any column of L,
221  * including the diagonal.
222  *
223  * Info [14..19] are not used in the current version, but may be used in
224  * future versions.
225  */
226 
227 /* ------------------------------------------------------------------------- */
228 /* direct interface to AMD */
229 /* ------------------------------------------------------------------------- */
230 
231 /* amesos_amd_2 is the primary AMD ordering routine. It is not meant to be
232  * user-callable because of its restrictive inputs and because it destroys
233  * the user's input matrix. It does not check its inputs for errors, either.
234  * However, if you can work with these restrictions it can be faster than
235  * amesos_amd_order and use less memory (assuming that you can create your own copy
236  * of the matrix for AMD to destroy). Refer to AMD/Source/amesos_amd_2.c for a
237  * description of each parameter. */
238 
239 void amesos_amd_2
240 (
241  int n,
242  int Pe [ ],
243  int Iw [ ],
244  int Len [ ],
245  int iwlen,
246  int pfree,
247  int Nv [ ],
248  int Next [ ],
249  int Last [ ],
250  int Head [ ],
251  int Elen [ ],
252  int Degree [ ],
253  int W [ ],
254  double Control [ ],
255  double Info [ ]
256 ) ;
257 
258 void amesos_amd_l2
259 (
260  UF_long n,
261  UF_long Pe [ ],
262  UF_long Iw [ ],
263  UF_long Len [ ],
264  UF_long iwlen,
265  UF_long pfree,
266  UF_long Nv [ ],
267  UF_long Next [ ],
268  UF_long Last [ ],
269  UF_long Head [ ],
270  UF_long Elen [ ],
271  UF_long Degree [ ],
272  UF_long W [ ],
273  double Control [ ],
274  double Info [ ]
275 ) ;
276 
277 /* ------------------------------------------------------------------------- */
278 /* amesos_amd_valid */
279 /* ------------------------------------------------------------------------- */
280 
281 /* Returns AMD_OK or AMD_OK_BUT_JUMBLED if the matrix is valid as input to
282  * amd_order; the latter is returned if the matrix has unsorted and/or
283  * duplicate row indices in one or more columns. Returns AMD_INVALID if the
284  * matrix cannot be passed to amesos_amd_order. For amesos_amd_order, the matrix must also
285  * be square. The first two arguments are the number of rows and the number
286  * of columns of the matrix. For its use in AMD, these must both equal n.
287  *
288  * NOTE: this routine returned TRUE/FALSE in v1.2 and earlier.
289  */
290 
292 (
293  int n_row, /* # of rows */
294  int n_col, /* # of columns */
295  const int Ap [ ], /* column pointers, of size n_col+1 */
296  const int Ai [ ] /* row indices, of size Ap [n_col] */
297 ) ;
298 
300 (
301  UF_long n_row,
302  UF_long n_col,
303  const UF_long Ap [ ],
304  const UF_long Ai [ ]
305 ) ;
306 
307 /* ------------------------------------------------------------------------- */
308 /* AMD memory manager and printf routines */
309 /* ------------------------------------------------------------------------- */
310 
311 /* The user can redefine these to change the malloc, free, and printf routines
312  * that AMD uses. */
313 
314 #ifndef EXTERN
315 #define EXTERN extern
316 #endif
317 
318 EXTERN void *(*amesos_amd_malloc) (size_t) ; /* pointer to malloc */
319 EXTERN void (*amesos_amd_free) (void *) ; /* pointer to free */
320 EXTERN void *(*amesos_amd_realloc) (void *, size_t) ; /* pointer to realloc */
321 EXTERN void *(*amesos_amd_calloc) (size_t, size_t) ; /* pointer to calloc */
322 EXTERN int (*amesos_amd_printf) (const char *, ...) ; /* pointer to printf */
323 
324 /* ------------------------------------------------------------------------- */
325 /* AMD Control and Info arrays */
326 /* ------------------------------------------------------------------------- */
327 
328 /* amesos_amd_defaults: sets the default control settings */
329 void amesos_amd_defaults (double Control [ ]) ;
330 void amesos_amd_l_defaults (double Control [ ]) ;
331 
332 /* amesos_amd_control: prints the control settings */
333 void amesos_amd_control (double Control [ ]) ;
334 void amesos_amd_l_control (double Control [ ]) ;
335 
336 /* amesos_amd_info: prints the statistics */
337 void amesos_amd_info (double Info [ ]) ;
338 void amesos_amd_l_info (double Info [ ]) ;
339 
340 #define AMD_CONTROL 5 /* size of Control array */
341 #define AMD_INFO 20 /* size of Info array */
342 
343 /* contents of Control */
344 #define AMD_DENSE 0 /* "dense" if degree > Control [0] * sqrt (n) */
345 #define AMD_AGGRESSIVE 1 /* do aggressive absorption if Control [1] != 0 */
346 
347 /* default Control settings */
348 #define AMD_DEFAULT_DENSE 10.0 /* default "dense" degree 10*sqrt(n) */
349 #define AMD_DEFAULT_AGGRESSIVE 1 /* do aggressive absorption by default */
350 
351 /* contents of Info */
352 #define AMD_STATUS 0 /* return value of amd_order and amd_l_order */
353 #define AMD_N 1 /* A is n-by-n */
354 #define AMD_NZ 2 /* number of nonzeros in A */
355 #define AMD_SYMMETRY 3 /* symmetry of pattern (1 is sym., 0 is unsym.) */
356 #define AMD_NZDIAG 4 /* # of entries on diagonal */
357 #define AMD_NZ_A_PLUS_AT 5 /* nz in A+A' */
358 #define AMD_NDENSE 6 /* number of "dense" rows/columns in A */
359 #define AMD_MEMORY 7 /* amount of memory used by AMD */
360 #define AMD_NCMPA 8 /* number of garbage collections in AMD */
361 #define AMD_LNZ 9 /* approx. nz in L, excluding the diagonal */
362 #define AMD_NDIV 10 /* number of fl. point divides for LU and LDL' */
363 #define AMD_NMULTSUBS_LDL 11 /* number of fl. point (*,-) pairs for LDL' */
364 #define AMD_NMULTSUBS_LU 12 /* number of fl. point (*,-) pairs for LU */
365 #define AMD_DMAX 13 /* max nz. in any column of L, incl. diagonal */
366 
367 /* ------------------------------------------------------------------------- */
368 /* return values of AMD */
369 /* ------------------------------------------------------------------------- */
370 
371 #define AMD_OK 0 /* success */
372 #define AMD_OUT_OF_MEMORY -1 /* malloc failed, or problem too large */
373 #define AMD_INVALID -2 /* input arguments are not valid */
374 #define AMD_OK_BUT_JUMBLED 1 /* input matrix is OK for amd_order, but
375  * columns were not sorted, and/or duplicate entries were present. AMD had
376  * to do extra work before ordering the matrix. This is a warning, not an
377  * error. */
378 
379 /* ========================================================================== */
380 /* === AMD version ========================================================== */
381 /* ========================================================================== */
382 
383 /* AMD Version 1.2 and later include the following definitions.
384  * As an example, to test if the version you are using is 1.2 or later:
385  *
386  * #ifdef AMD_VERSION
387  * if (AMD_VERSION >= AMD_VERSION_CODE (1,2)) ...
388  * #endif
389  *
390  * This also works during compile-time:
391  *
392  * #if defined(AMD_VERSION) && (AMD_VERSION >= AMD_VERSION_CODE (1,2))
393  * printf ("This is version 1.2 or later\n") ;
394  * #else
395  * printf ("This is an early version\n") ;
396  * #endif
397  *
398  * Versions 1.1 and earlier of AMD do not include a #define'd version number.
399  */
400 
401 #define AMD_DATE "May 31, 2007"
402 #define AMD_VERSION_CODE(main,sub) ((main) * 1000 + (sub))
403 #define AMD_MAIN_VERSION 2
404 #define AMD_SUB_VERSION 2
405 #define AMD_SUBSUB_VERSION 0
406 #define AMD_VERSION AMD_VERSION_CODE(AMD_MAIN_VERSION,AMD_SUB_VERSION)
407 
408 #ifdef __cplusplus
409 }
410 #endif
411 
412 #endif
UF_long amesos_amd_l_order(UF_long n, const UF_long Ap[], const UF_long Ai[], UF_long P[], double Control[], double Info[])
void amesos_amd_l_info(double Info[])
void amesos_amd_l_defaults(double Control[])
void amesos_amd_control(double Control[])
#define P(k)
EXTERN void(* amesos_amd_free)(void *)
Definition: amesos_amd.h:319
#define EXTERN
Definition: amesos_amd.h:315
int amesos_amd_order(int n, const int Ap[], const int Ai[], int P[], double Control[], double Info[])
EXTERN int(* amesos_amd_printf)(const char *,...)
Definition: amesos_amd.h:322
void amesos_amd_info(double Info[])
int amesos_amd_valid(int n_row, int n_col, const int Ap[], const int Ai[])
void amesos_amd_defaults(double Control[])
UF_long amesos_amd_l_valid(UF_long n_row, UF_long n_col, const UF_long Ap[], const UF_long Ai[])
void amesos_amd_2(int n, int Pe[], int Iw[], int Len[], int iwlen, int pfree, int Nv[], int Next[], int Last[], int Head[], int Elen[], int Degree[], int W[], double Control[], double Info[])
void amesos_amd_l2(UF_long n, UF_long Pe[], UF_long Iw[], UF_long Len[], UF_long iwlen, UF_long pfree, UF_long Nv[], UF_long Next[], UF_long Last[], UF_long Head[], UF_long Elen[], UF_long Degree[], UF_long W[], double Control[], double Info[])
void amesos_amd_l_control(double Control[])
#define UF_long
int n