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