Amesos Package Browser (Single Doxygen Collection)  Development
 All Classes Namespaces Files Functions Variables Typedefs Enumerations Enumerator Friends Macros Groups Pages
amesos_cholmod_cholesky.h
Go to the documentation of this file.
1 /* ========================================================================== */
2 /* === Include/cholmod_cholesky.h =========================================== */
3 /* ========================================================================== */
4 
5 /* -----------------------------------------------------------------------------
6  * CHOLMOD/Include/cholmod_cholesky.h. Copyright (C) 2005-2006, Timothy A. Davis
7  * CHOLMOD/Include/cholmod_cholesky.h is licensed under Version 2.1 of the GNU
8  * Lesser General Public License. See lesser.txt for a text of the license.
9  * CHOLMOD is also available under other licenses; contact authors for details.
10  * http://www.cise.ufl.edu/research/sparse
11  * -------------------------------------------------------------------------- */
12 
13 /* CHOLMOD Cholesky module.
14  *
15  * Sparse Cholesky routines: analysis, factorization, and solve.
16  *
17  * The primary routines are all that a user requires to order, analyze, and
18  * factorize a sparse symmetric positive definite matrix A (or A*A'), and
19  * to solve Ax=b (or A*A'x=b). The primary routines rely on the secondary
20  * routines, the CHOLMOD Core module, and the AMD and COLAMD packages. They
21  * make optional use of the CHOLMOD Supernodal and Partition modules, the
22  * METIS package, and the CCOLAMD package.
23  *
24  * Primary routines:
25  * -----------------
26  *
27  * cholmod_analyze order and analyze (simplicial or supernodal)
28  * cholmod_factorize simplicial or supernodal Cholesky factorization
29  * cholmod_solve solve a linear system (simplicial or supernodal)
30  * cholmod_spsolve solve a linear system (sparse x and b)
31  *
32  * Secondary routines:
33  * ------------------
34  *
35  * cholmod_analyze_p analyze, with user-provided permutation or f set
36  * cholmod_factorize_p factorize, with user-provided permutation or f
37  * cholmod_analyze_ordering analyze a fill-reducing ordering
38  * cholmod_etree find the elimination tree
39  * cholmod_rowcolcounts compute the row/column counts of L
40  * cholmod_amd order using AMD
41  * cholmod_colamd order using COLAMD
42  * cholmod_rowfac incremental simplicial factorization
43  * cholmod_rowfac_mask rowfac, specific to LPDASA
44  * cholmod_row_subtree find the nonzero pattern of a row of L
45  * cholmod_resymbol recompute the symbolic pattern of L
46  * cholmod_resymbol_noperm recompute the symbolic pattern of L, no L->Perm
47  * cholmod_postorder postorder a tree
48  *
49  * Requires the Core module, and two packages: AMD and COLAMD.
50  * Optionally uses the Supernodal and Partition modules.
51  * Required by the Partition module.
52  */
53 
54 #ifndef AMESOS_CHOLMOD_CHOLESKY_H
55 #define AMESOS_CHOLMOD_CHOLESKY_H
56 
57 #include "amesos_cholmod_config.h"
58 #include "amesos_cholmod_core.h"
59 
60 #ifndef NPARTITION
62 #endif
63 
64 /* -------------------------------------------------------------------------- */
65 /* cholmod_analyze: order and analyze (simplicial or supernodal) */
66 /* -------------------------------------------------------------------------- */
67 
68 /* Orders and analyzes A, AA', PAP', or PAA'P' and returns a symbolic factor
69  * that can later be passed to cholmod_factorize. */
70 
72 (
73  /* ---- input ---- */
74  cholmod_sparse *A, /* matrix to order and analyze */
75  /* --------------- */
76  cholmod_common *Common
77 ) ;
78 
80 
81 /* -------------------------------------------------------------------------- */
82 /* cholmod_analyze_p: analyze, with user-provided permutation or f set */
83 /* -------------------------------------------------------------------------- */
84 
85 /* Orders and analyzes A, AA', PAP', PAA'P', FF', or PFF'P and returns a
86  * symbolic factor that can later be passed to cholmod_factorize, where
87  * F = A(:,fset) if fset is not NULL and A->stype is zero.
88  * UserPerm is tried if non-NULL. */
89 
91 (
92  /* ---- input ---- */
93  cholmod_sparse *A, /* matrix to order and analyze */
94  int *UserPerm, /* user-provided permutation, size A->nrow */
95  int *fset, /* subset of 0:(A->ncol)-1 */
96  size_t fsize, /* size of fset */
97  /* --------------- */
98  cholmod_common *Common
99 ) ;
100 
102  size_t, cholmod_common *) ;
103 
104 /* -------------------------------------------------------------------------- */
105 /* cholmod_factorize: simplicial or supernodal Cholesky factorization */
106 /* -------------------------------------------------------------------------- */
107 
108 /* Factorizes PAP' (or PAA'P' if A->stype is 0), using a factor obtained
109  * from cholmod_analyze. The analysis can be re-used simply by calling this
110  * routine a second time with another matrix. A must have the same nonzero
111  * pattern as that passed to cholmod_analyze. */
112 
114 (
115  /* ---- input ---- */
116  cholmod_sparse *A, /* matrix to factorize */
117  /* ---- in/out --- */
118  cholmod_factor *L, /* resulting factorization */
119  /* --------------- */
120  cholmod_common *Common
121 ) ;
122 
124 
125 /* -------------------------------------------------------------------------- */
126 /* cholmod_factorize_p: factorize, with user-provided permutation or fset */
127 /* -------------------------------------------------------------------------- */
128 
129 /* Same as cholmod_factorize, but with more options. */
130 
132 (
133  /* ---- input ---- */
134  cholmod_sparse *A, /* matrix to factorize */
135  double beta [2], /* factorize beta*I+A or beta*I+A'*A */
136  int *fset, /* subset of 0:(A->ncol)-1 */
137  size_t fsize, /* size of fset */
138  /* ---- in/out --- */
139  cholmod_factor *L, /* resulting factorization */
140  /* --------------- */
141  cholmod_common *Common
142 ) ;
143 
144 int amesos_cholmod_l_factorize_p (cholmod_sparse *, double *, UF_long *, size_t,
146 
147 /* -------------------------------------------------------------------------- */
148 /* cholmod_solve: solve a linear system (simplicial or supernodal) */
149 /* -------------------------------------------------------------------------- */
150 
151 /* Solves one of many linear systems with a dense right-hand-side, using the
152  * factorization from cholmod_factorize (or as modified by any other CHOLMOD
153  * routine). D is identity for LL' factorizations. */
154 
155 #define CHOLMOD_A 0 /* solve Ax=b */
156 #define CHOLMOD_LDLt 1 /* solve LDL'x=b */
157 #define CHOLMOD_LD 2 /* solve LDx=b */
158 #define CHOLMOD_DLt 3 /* solve DL'x=b */
159 #define CHOLMOD_L 4 /* solve Lx=b */
160 #define CHOLMOD_Lt 5 /* solve L'x=b */
161 #define CHOLMOD_D 6 /* solve Dx=b */
162 #define CHOLMOD_P 7 /* permute x=Px */
163 #define CHOLMOD_Pt 8 /* permute x=P'x */
164 
165 cholmod_dense *amesos_cholmod_solve /* returns the solution X */
166 (
167  /* ---- input ---- */
168  int sys, /* system to solve */
169  cholmod_factor *L, /* factorization to use */
170  cholmod_dense *B, /* right-hand-side */
171  /* --------------- */
172  cholmod_common *Common
173 ) ;
174 
176  cholmod_common *) ;
177 
178 /* -------------------------------------------------------------------------- */
179 /* cholmod_spsolve: solve a linear system with a sparse right-hand-side */
180 /* -------------------------------------------------------------------------- */
181 
183 (
184  /* ---- input ---- */
185  int sys, /* system to solve */
186  cholmod_factor *L, /* factorization to use */
187  cholmod_sparse *B, /* right-hand-side */
188  /* --------------- */
189  cholmod_common *Common
190 ) ;
191 
193  cholmod_common *) ;
194 
195 /* -------------------------------------------------------------------------- */
196 /* cholmod_etree: find the elimination tree of A or A'*A */
197 /* -------------------------------------------------------------------------- */
198 
200 (
201  /* ---- input ---- */
202  cholmod_sparse *A,
203  /* ---- output --- */
204  int *Parent, /* size ncol. Parent [j] = p if p is the parent of j */
205  /* --------------- */
206  cholmod_common *Common
207 ) ;
208 
210 
211 /* -------------------------------------------------------------------------- */
212 /* cholmod_rowcolcounts: compute the row/column counts of L */
213 /* -------------------------------------------------------------------------- */
214 
216 (
217  /* ---- input ---- */
218  cholmod_sparse *A, /* matrix to analyze */
219  int *fset, /* subset of 0:(A->ncol)-1 */
220  size_t fsize, /* size of fset */
221  int *Parent, /* size nrow. Parent [i] = p if p is the parent of i */
222  int *Post, /* size nrow. Post [k] = i if i is the kth node in
223  * the postordered etree. */
224  /* ---- output --- */
225  int *RowCount, /* size nrow. RowCount [i] = # entries in the ith row of
226  * L, including the diagonal. */
227  int *ColCount, /* size nrow. ColCount [i] = # entries in the ith
228  * column of L, including the diagonal. */
229  int *First, /* size nrow. First [i] = k is the least postordering
230  * of any descendant of i. */
231  int *Level, /* size nrow. Level [i] is the length of the path from
232  * i to the root, with Level [root] = 0. */
233  /* --------------- */
234  cholmod_common *Common
235 ) ;
236 
239 
240 /* -------------------------------------------------------------------------- */
241 /* cholmod_analyze_ordering: analyze a fill-reducing ordering */
242 /* -------------------------------------------------------------------------- */
243 
245 (
246  /* ---- input ---- */
247  cholmod_sparse *A, /* matrix to analyze */
248  int ordering, /* ordering method used */
249  int *Perm, /* size n, fill-reducing permutation to analyze */
250  int *fset, /* subset of 0:(A->ncol)-1 */
251  size_t fsize, /* size of fset */
252  /* ---- output --- */
253  int *Parent, /* size n, elimination tree */
254  int *Post, /* size n, postordering of elimination tree */
255  int *ColCount, /* size n, nnz in each column of L */
256  /* ---- workspace */
257  int *First, /* size nworkspace for cholmod_postorder */
258  int *Level, /* size n workspace for cholmod_postorder */
259  /* --------------- */
260  cholmod_common *Common
261 ) ;
262 
264  size_t, UF_long *, UF_long *, UF_long *, UF_long *, UF_long *,
265  cholmod_common *) ;
266 
267 /* -------------------------------------------------------------------------- */
268 /* cholmod_amd: order using AMD */
269 /* -------------------------------------------------------------------------- */
270 
271 /* Finds a permutation P to reduce fill-in in the factorization of P*A*P'
272  * or P*A*A'P' */
273 
275 (
276  /* ---- input ---- */
277  cholmod_sparse *A, /* matrix to order */
278  int *fset, /* subset of 0:(A->ncol)-1 */
279  size_t fsize, /* size of fset */
280  /* ---- output --- */
281  int *Perm, /* size A->nrow, output permutation */
282  /* --------------- */
283  cholmod_common *Common
284 ) ;
285 
287  cholmod_common *) ;
288 
289 /* -------------------------------------------------------------------------- */
290 /* cholmod_colamd: order using COLAMD */
291 /* -------------------------------------------------------------------------- */
292 
293 /* Finds a permutation P to reduce fill-in in the factorization of P*A*A'*P'.
294  * Orders F*F' where F = A (:,fset) if fset is not NULL */
295 
297 (
298  /* ---- input ---- */
299  cholmod_sparse *A, /* matrix to order */
300  int *fset, /* subset of 0:(A->ncol)-1 */
301  size_t fsize, /* size of fset */
302  int postorder, /* if TRUE, follow with a coletree postorder */
303  /* ---- output --- */
304  int *Perm, /* size A->nrow, output permutation */
305  /* --------------- */
306  cholmod_common *Common
307 ) ;
308 
309 int amesos_cholmod_l_colamd (cholmod_sparse *, UF_long *, size_t, int, UF_long *,
310  cholmod_common *) ;
311 
312 /* -------------------------------------------------------------------------- */
313 /* cholmod_rowfac: incremental simplicial factorization */
314 /* -------------------------------------------------------------------------- */
315 
316 /* Partial or complete simplicial factorization. Rows and columns kstart:kend-1
317  * of L and D must be initially equal to rows/columns kstart:kend-1 of the
318  * identity matrix. Row k can only be factorized if all descendants of node
319  * k in the elimination tree have been factorized. */
320 
322 (
323  /* ---- input ---- */
324  cholmod_sparse *A, /* matrix to factorize */
325  cholmod_sparse *F, /* used for A*A' case only. F=A' or A(:,fset)' */
326  double beta [2], /* factorize beta*I+A or beta*I+A'*A */
327  size_t kstart, /* first row to factorize */
328  size_t kend, /* last row to factorize is kend-1 */
329  /* ---- in/out --- */
330  cholmod_factor *L,
331  /* --------------- */
332  cholmod_common *Common
333 ) ;
334 
335 int amesos_cholmod_l_rowfac (cholmod_sparse *, cholmod_sparse *, double *, size_t,
336  size_t, cholmod_factor *, cholmod_common *) ;
337 
338 /* -------------------------------------------------------------------------- */
339 /* cholmod_rowfac_mask: incremental simplicial factorization */
340 /* -------------------------------------------------------------------------- */
341 
342 /* cholmod_rowfac_mask is a version of cholmod_rowfac that is specific to
343  * LPDASA. It is unlikely to be needed by any other application. */
344 
346 (
347  /* ---- input ---- */
348  cholmod_sparse *A, /* matrix to factorize */
349  cholmod_sparse *F, /* used for A*A' case only. F=A' or A(:,fset)' */
350  double beta [2], /* factorize beta*I+A or beta*I+A'*A */
351  size_t kstart, /* first row to factorize */
352  size_t kend, /* last row to factorize is kend-1 */
353  int *mask, /* if mask[i] >= 0, then set row i to zero */
354  int *RLinkUp, /* link list of rows to compute */
355  /* ---- in/out --- */
356  cholmod_factor *L,
357  /* --------------- */
358  cholmod_common *Common
359 ) ;
360 
362  size_t, UF_long *, UF_long *, cholmod_factor *, cholmod_common *) ;
363 
364 /* -------------------------------------------------------------------------- */
365 /* cholmod_row_subtree: find the nonzero pattern of a row of L */
366 /* -------------------------------------------------------------------------- */
367 
368 /* Find the nonzero pattern of x for the system Lx=b where L = (0:k-1,0:k-1)
369  * and b = kth column of A or A*A' (rows 0 to k-1 only) */
370 
372 (
373  /* ---- input ---- */
374  cholmod_sparse *A, /* matrix to analyze */
375  cholmod_sparse *F, /* used for A*A' case only. F=A' or A(:,fset)' */
376  size_t k, /* row k of L */
377  int *Parent, /* elimination tree */
378  /* ---- output --- */
379  cholmod_sparse *R, /* pattern of L(k,:), 1-by-n with R->nzmax >= n */
380  /* --------------- */
381  cholmod_common *Common
382 ) ;
383 
386 
387 /* -------------------------------------------------------------------------- */
388 /* cholmod_row_lsubtree: find the nonzero pattern of a row of L */
389 /* -------------------------------------------------------------------------- */
390 
391 /* Identical to cholmod_row_subtree, except that it finds the elimination tree
392  * from L itself. */
393 
395 (
396  /* ---- input ---- */
397  cholmod_sparse *A, /* matrix to analyze */
398  int *Fi, size_t fnz, /* nonzero pattern of kth row of A', not required
399  * for the symmetric case. Need not be sorted. */
400  size_t k, /* row k of L */
401  cholmod_factor *L, /* the factor L from which parent(i) is derived */
402  /* ---- output --- */
403  cholmod_sparse *R, /* pattern of L(k,:), 1-by-n with R->nzmax >= n */
404  /* --------------- */
405  cholmod_common *Common
406 ) ;
407 
409  size_t, cholmod_factor *, cholmod_sparse *, cholmod_common *) ;
410 
411 /* -------------------------------------------------------------------------- */
412 /* cholmod_resymbol: recompute the symbolic pattern of L */
413 /* -------------------------------------------------------------------------- */
414 
415 /* Remove entries from L that are not in the factorization of P*A*P', P*A*A'*P',
416  * or P*F*F'*P' (depending on A->stype and whether fset is NULL or not).
417  *
418  * cholmod_resymbol is the same as cholmod_resymbol_noperm, except that it
419  * first permutes A according to L->Perm. A can be upper/lower/unsymmetric,
420  * in contrast to cholmod_resymbol_noperm (which can be lower or unsym). */
421 
423 (
424  /* ---- input ---- */
425  cholmod_sparse *A, /* matrix to analyze */
426  int *fset, /* subset of 0:(A->ncol)-1 */
427  size_t fsize, /* size of fset */
428  int pack, /* if TRUE, pack the columns of L */
429  /* ---- in/out --- */
430  cholmod_factor *L, /* factorization, entries pruned on output */
431  /* --------------- */
432  cholmod_common *Common
433 ) ;
434 
435 int amesos_cholmod_l_resymbol (cholmod_sparse *, UF_long *, size_t, int,
437 
438 /* -------------------------------------------------------------------------- */
439 /* cholmod_resymbol_noperm: recompute the symbolic pattern of L, no L->Perm */
440 /* -------------------------------------------------------------------------- */
441 
442 /* Remove entries from L that are not in the factorization of A, A*A',
443  * or F*F' (depending on A->stype and whether fset is NULL or not). */
444 
446 (
447  /* ---- input ---- */
448  cholmod_sparse *A, /* matrix to analyze */
449  int *fset, /* subset of 0:(A->ncol)-1 */
450  size_t fsize, /* size of fset */
451  int pack, /* if TRUE, pack the columns of L */
452  /* ---- in/out --- */
453  cholmod_factor *L, /* factorization, entries pruned on output */
454  /* --------------- */
455  cholmod_common *Common
456 ) ;
457 
460 
461 /* -------------------------------------------------------------------------- */
462 /* cholmod_rcond: compute rough estimate of reciprocal of condition number */
463 /* -------------------------------------------------------------------------- */
464 
465 double amesos_cholmod_rcond /* return min(diag(L)) / max(diag(L)) */
466 (
467  /* ---- input ---- */
468  cholmod_factor *L,
469  /* --------------- */
470  cholmod_common *Common
471 ) ;
472 
474 
475 /* -------------------------------------------------------------------------- */
476 /* cholmod_postorder: Compute the postorder of a tree */
477 /* -------------------------------------------------------------------------- */
478 
479 UF_long amesos_cholmod_postorder /* return # of nodes postordered */
480 (
481  /* ---- input ---- */
482  int *Parent, /* size n. Parent [j] = p if p is the parent of j */
483  size_t n,
484  int *Weight_p, /* size n, optional. Weight [j] is weight of node j */
485  /* ---- output --- */
486  int *Post, /* size n. Post [k] = j is kth in postordered tree */
487  /* --------------- */
488  cholmod_common *Common
489 ) ;
490 
492  cholmod_common *) ;
493 
494 #endif
cholmod_sparse * amesos_cholmod_spsolve(int sys, cholmod_factor *L, cholmod_sparse *B, cholmod_common *Common)
int amesos_cholmod_l_rowfac(cholmod_sparse *, cholmod_sparse *, double *, size_t, size_t, cholmod_factor *, cholmod_common *)
int amesos_cholmod_resymbol(cholmod_sparse *A, int *fset, size_t fsize, int pack, cholmod_factor *L, cholmod_common *Common)
int amesos_cholmod_l_row_subtree(cholmod_sparse *, cholmod_sparse *, size_t, UF_long *, cholmod_sparse *, cholmod_common *)
int amesos_cholmod_l_rowcolcounts(cholmod_sparse *, UF_long *, size_t, UF_long *, UF_long *, UF_long *, UF_long *, UF_long *, UF_long *, cholmod_common *)
int amesos_cholmod_row_subtree(cholmod_sparse *A, cholmod_sparse *F, size_t k, int *Parent, cholmod_sparse *R, cholmod_common *Common)
int amesos_cholmod_rowfac(cholmod_sparse *A, cholmod_sparse *F, double beta[2], size_t kstart, size_t kend, cholmod_factor *L, cholmod_common *Common)
int amesos_cholmod_row_lsubtree(cholmod_sparse *A, int *Fi, size_t fnz, size_t k, cholmod_factor *L, cholmod_sparse *R, cholmod_common *Common)
int amesos_cholmod_l_etree(cholmod_sparse *, UF_long *, cholmod_common *)
int amesos_cholmod_analyze_ordering(cholmod_sparse *A, int ordering, int *Perm, int *fset, size_t fsize, int *Parent, int *Post, int *ColCount, int *First, int *Level, cholmod_common *Common)
int amesos_cholmod_l_resymbol(cholmod_sparse *, UF_long *, size_t, int, cholmod_factor *, cholmod_common *)
int amesos_cholmod_l_rowfac_mask(cholmod_sparse *, cholmod_sparse *, double *, size_t, size_t, UF_long *, UF_long *, cholmod_factor *, cholmod_common *)
cholmod_factor * cholmod_l_analyze(cholmod_sparse *, cholmod_common *)
cholmod_factor * amesos_cholmod_analyze(cholmod_sparse *A, cholmod_common *Common)
int amesos_cholmod_l_resymbol_noperm(cholmod_sparse *, UF_long *, size_t, int, cholmod_factor *, cholmod_common *)
cholmod_dense * amesos_cholmod_l_solve(int, cholmod_factor *, cholmod_dense *, cholmod_common *)
cholmod_dense * amesos_cholmod_solve(int sys, cholmod_factor *L, cholmod_dense *B, cholmod_common *Common)
int amesos_cholmod_l_factorize(cholmod_sparse *, cholmod_factor *, cholmod_common *)
int amesos_cholmod_rowcolcounts(cholmod_sparse *A, int *fset, size_t fsize, int *Parent, int *Post, int *RowCount, int *ColCount, int *First, int *Level, cholmod_common *Common)
int amesos_cholmod_colamd(cholmod_sparse *A, int *fset, size_t fsize, int postorder, int *Perm, cholmod_common *Common)
int amesos_cholmod_l_analyze_ordering(cholmod_sparse *, int, UF_long *, UF_long *, size_t, UF_long *, UF_long *, UF_long *, UF_long *, UF_long *, cholmod_common *)
int amesos_cholmod_l_factorize_p(cholmod_sparse *, double *, UF_long *, size_t, cholmod_factor *, cholmod_common *)
UF_long amesos_cholmod_l_postorder(UF_long *, size_t, UF_long *, UF_long *, cholmod_common *)
int amesos_cholmod_rowfac_mask(cholmod_sparse *A, cholmod_sparse *F, double beta[2], size_t kstart, size_t kend, int *mask, int *RLinkUp, cholmod_factor *L, cholmod_common *Common)
int amesos_cholmod_l_row_lsubtree(cholmod_sparse *, UF_long *, size_t, size_t, cholmod_factor *, cholmod_sparse *, cholmod_common *)
double amesos_cholmod_l_rcond(cholmod_factor *, cholmod_common *)
int amesos_cholmod_factorize_p(cholmod_sparse *A, double beta[2], int *fset, size_t fsize, cholmod_factor *L, cholmod_common *Common)
int amesos_cholmod_etree(cholmod_sparse *A, int *Parent, cholmod_common *Common)
UF_long amesos_cholmod_postorder(int *Parent, size_t n, int *Weight_p, int *Post, cholmod_common *Common)
cholmod_factor * amesos_cholmod_l_analyze_p(cholmod_sparse *, UF_long *, UF_long *, size_t, cholmod_common *)
int amesos_cholmod_resymbol_noperm(cholmod_sparse *A, int *fset, size_t fsize, int pack, cholmod_factor *L, cholmod_common *Common)
int amesos_cholmod_l_amd(cholmod_sparse *, UF_long *, size_t, UF_long *, cholmod_common *)
int amesos_cholmod_l_colamd(cholmod_sparse *, UF_long *, size_t, int, UF_long *, cholmod_common *)
int amesos_cholmod_amd(cholmod_sparse *A, int *fset, size_t fsize, int *Perm, cholmod_common *Common)
double amesos_cholmod_rcond(cholmod_factor *L, cholmod_common *Common)
int amesos_cholmod_factorize(cholmod_sparse *A, cholmod_factor *L, cholmod_common *Common)
#define UF_long
cholmod_sparse * amesos_cholmod_l_spsolve(int, cholmod_factor *, cholmod_sparse *, cholmod_common *)
int n
cholmod_factor * amesos_cholmod_analyze_p(cholmod_sparse *A, int *UserPerm, int *fset, size_t fsize, cholmod_common *Common)