Actual source code: sbaijcholmod.c
1: /*
2: Provides an interface to the CHOLMOD sparse solver available through SuiteSparse version 4.2.1
4: When built with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the
5: integer type in UMFPACK, otherwise it will use int. This means
6: all integers in this file as simply declared as PetscInt. Also it means
7: that one cannot use 64BIT_INDICES on 32-bit pointer systems [as Suitesparse_long is 32-bit only]
9: */
11: #include <../src/mat/impls/sbaij/seq/sbaij.h>
12: #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h>
14: /*
15: This is a terrible hack, but it allows the error handler to retain a context.
16: Note that this hack really cannot be made both reentrant and concurrent.
17: */
18: static Mat static_F;
20: static void CholmodErrorHandler(int status, const char *file, int line, const char *message)
21: {
22: PetscFunctionBegin;
23: if (status > CHOLMOD_OK) {
24: PetscCallVoid(PetscInfo(static_F, "CHOLMOD warning %d at %s:%d: %s\n", status, file, line, message));
25: } else if (status == CHOLMOD_OK) { /* Documentation says this can happen, but why? */
26: PetscCallVoid(PetscInfo(static_F, "CHOLMOD OK at %s:%d: %s\n", file, line, message));
27: } else {
28: PetscCallVoid(PetscErrorPrintf("CHOLMOD error %d at %s:%d: %s\n", status, file, line, message));
29: }
30: PetscFunctionReturnVoid();
31: }
33: #define CHOLMOD_OPTION_DOUBLE(name, help) \
34: do { \
35: PetscReal tmp = (PetscReal)c->name; \
36: PetscCall(PetscOptionsReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \
37: c->name = (double)tmp; \
38: } while (0)
40: #define CHOLMOD_OPTION_INT(name, help) \
41: do { \
42: PetscInt tmp = (PetscInt)c->name; \
43: PetscCall(PetscOptionsInt("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \
44: c->name = (int)tmp; \
45: } while (0)
47: #define CHOLMOD_OPTION_SIZE_T(name, help) \
48: do { \
49: PetscReal tmp = (PetscInt)c->name; \
50: PetscCall(PetscOptionsBoundedReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL, 0.0)); \
51: c->name = (size_t)tmp; \
52: } while (0)
54: #define CHOLMOD_OPTION_BOOL(name, help) \
55: do { \
56: PetscBool tmp = (PetscBool)!!c->name; \
57: PetscCall(PetscOptionsBool("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \
58: c->name = (int)tmp; \
59: } while (0)
61: static PetscErrorCode CholmodSetOptions(Mat F)
62: {
63: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
64: cholmod_common *c = chol->common;
65: PetscBool flg;
67: PetscFunctionBegin;
68: PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "CHOLMOD Options", "Mat");
69: CHOLMOD_OPTION_INT(nmethods, "Number of different ordering methods to try");
71: #if defined(PETSC_USE_SUITESPARSE_GPU)
72: c->useGPU = 1;
73: CHOLMOD_OPTION_INT(useGPU, "Use GPU for BLAS 1, otherwise 0");
74: CHOLMOD_OPTION_SIZE_T(maxGpuMemBytes, "Maximum memory to allocate on the GPU");
75: CHOLMOD_OPTION_DOUBLE(maxGpuMemFraction, "Fraction of available GPU memory to allocate");
76: #endif
78: /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */
79: chol->pack = (PetscBool)c->final_pack;
80: PetscCall(PetscOptionsBool("-mat_cholmod_pack", "Pack factors after factorization [disable for frequent repeat factorization]", "None", chol->pack, &chol->pack, NULL));
81: c->final_pack = (int)chol->pack;
83: CHOLMOD_OPTION_DOUBLE(dbound, "Minimum absolute value of diagonal entries of D");
84: CHOLMOD_OPTION_DOUBLE(grow0, "Global growth ratio when factors are modified");
85: CHOLMOD_OPTION_DOUBLE(grow1, "Column growth ratio when factors are modified");
86: CHOLMOD_OPTION_SIZE_T(grow2, "Affine column growth constant when factors are modified");
87: CHOLMOD_OPTION_SIZE_T(maxrank, "Max rank of update, larger values are faster but use more memory [2,4,8]");
88: {
89: static const char *const list[] = {"SIMPLICIAL", "AUTO", "SUPERNODAL", "MatCholmodFactorType", "MAT_CHOLMOD_FACTOR_", 0};
90: PetscCall(PetscOptionsEnum("-mat_cholmod_factor", "Factorization method", "None", list, (PetscEnum)c->supernodal, (PetscEnum *)&c->supernodal, NULL));
91: }
92: if (c->supernodal) CHOLMOD_OPTION_DOUBLE(supernodal_switch, "flop/nnz_L threshold for switching to supernodal factorization");
93: CHOLMOD_OPTION_BOOL(final_asis, "Leave factors \"as is\"");
94: CHOLMOD_OPTION_BOOL(final_pack, "Pack the columns when finished (use FALSE if the factors will be updated later)");
95: if (!c->final_asis) {
96: CHOLMOD_OPTION_BOOL(final_super, "Leave supernodal factors instead of converting to simplicial");
97: CHOLMOD_OPTION_BOOL(final_ll, "Turn LDL' factorization into LL'");
98: CHOLMOD_OPTION_BOOL(final_monotonic, "Ensure columns are monotonic when done");
99: CHOLMOD_OPTION_BOOL(final_resymbol, "Remove numerically zero values resulting from relaxed supernodal amalgamation");
100: }
101: {
102: PetscReal tmp[] = {(PetscReal)c->zrelax[0], (PetscReal)c->zrelax[1], (PetscReal)c->zrelax[2]};
103: PetscInt n = 3;
104: PetscCall(PetscOptionsRealArray("-mat_cholmod_zrelax", "3 real supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg));
105: PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_zrelax");
106: if (flg)
107: while (n--) c->zrelax[n] = (double)tmp[n];
108: }
109: {
110: PetscInt n, tmp[] = {(PetscInt)c->nrelax[0], (PetscInt)c->nrelax[1], (PetscInt)c->nrelax[2]};
111: PetscCall(PetscOptionsIntArray("-mat_cholmod_nrelax", "3 size_t supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg));
112: PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_nrelax");
113: if (flg)
114: while (n--) c->nrelax[n] = (size_t)tmp[n];
115: }
116: CHOLMOD_OPTION_BOOL(prefer_upper, "Work with upper triangular form [faster when using fill-reducing ordering, slower in natural ordering]");
117: CHOLMOD_OPTION_BOOL(default_nesdis, "Use NESDIS instead of METIS for nested dissection");
118: CHOLMOD_OPTION_INT(print, "Verbosity level");
119: PetscOptionsEnd();
120: PetscFunctionReturn(PETSC_SUCCESS);
121: }
123: PetscErrorCode CholmodStart(Mat F)
124: {
125: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
126: cholmod_common *c;
128: PetscFunctionBegin;
129: if (chol->common) PetscFunctionReturn(PETSC_SUCCESS);
130: PetscCall(PetscMalloc1(1, &chol->common));
131: PetscCallExternal(!cholmod_X_start, chol->common);
133: c = chol->common;
134: c->error_handler = CholmodErrorHandler;
135: PetscFunctionReturn(PETSC_SUCCESS);
136: }
138: static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc)
139: {
140: Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ *)A->data;
141: PetscBool vallocin = PETSC_FALSE;
143: PetscFunctionBegin;
144: PetscCall(PetscMemzero(C, sizeof(*C)));
145: /* CHOLMOD uses column alignment, SBAIJ stores the upper factor, so we pass it on as a lower factor, swapping the meaning of row and column */
146: C->nrow = (size_t)A->cmap->n;
147: C->ncol = (size_t)A->rmap->n;
148: C->nzmax = (size_t)sbaij->maxnz;
149: C->p = sbaij->i;
150: C->i = sbaij->j;
151: if (values) {
152: #if defined(PETSC_USE_COMPLEX)
153: /* we need to pass CHOLMOD the conjugate matrix */
154: PetscScalar *v;
156: PetscCall(PetscMalloc1(sbaij->maxnz, &v));
157: for (PetscInt i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]);
158: C->x = v;
159: vallocin = PETSC_TRUE;
160: #else
161: C->x = sbaij->a;
162: #endif
163: }
164: C->stype = -1;
165: C->itype = CHOLMOD_INT_TYPE;
166: C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
167: C->dtype = CHOLMOD_DOUBLE;
168: C->sorted = 1;
169: C->packed = 1;
170: *aijalloc = PETSC_FALSE;
171: *valloc = vallocin;
172: PetscFunctionReturn(PETSC_SUCCESS);
173: }
175: #define GET_ARRAY_READ 0
176: #define GET_ARRAY_WRITE 1
178: PetscErrorCode VecWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y)
179: {
180: PetscScalar *x;
181: PetscInt n;
183: PetscFunctionBegin;
184: PetscCall(PetscMemzero(Y, sizeof(*Y)));
185: switch (rw) {
186: case GET_ARRAY_READ:
187: PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
188: break;
189: case GET_ARRAY_WRITE:
190: PetscCall(VecGetArrayWrite(X, &x));
191: break;
192: default:
193: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
194: break;
195: }
196: PetscCall(VecGetSize(X, &n));
198: Y->x = x;
199: Y->nrow = n;
200: Y->ncol = 1;
201: Y->nzmax = n;
202: Y->d = n;
203: Y->xtype = CHOLMOD_SCALAR_TYPE;
204: Y->dtype = CHOLMOD_DOUBLE;
205: PetscFunctionReturn(PETSC_SUCCESS);
206: }
208: PetscErrorCode VecUnWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y)
209: {
210: PetscFunctionBegin;
211: switch (rw) {
212: case GET_ARRAY_READ:
213: PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&Y->x));
214: break;
215: case GET_ARRAY_WRITE:
216: PetscCall(VecRestoreArrayWrite(X, (PetscScalar **)&Y->x));
217: break;
218: default:
219: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
220: break;
221: }
222: PetscFunctionReturn(PETSC_SUCCESS);
223: }
225: PetscErrorCode MatDenseWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y)
226: {
227: PetscScalar *x;
228: PetscInt m, n, lda;
230: PetscFunctionBegin;
231: PetscCall(PetscMemzero(Y, sizeof(*Y)));
232: switch (rw) {
233: case GET_ARRAY_READ:
234: PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&x));
235: break;
236: case GET_ARRAY_WRITE:
237: PetscCall(MatDenseGetArrayWrite(X, &x));
238: break;
239: default:
240: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
241: break;
242: }
243: PetscCall(MatDenseGetLDA(X, &lda));
244: PetscCall(MatGetLocalSize(X, &m, &n));
246: Y->x = x;
247: Y->nrow = m;
248: Y->ncol = n;
249: Y->nzmax = lda * n;
250: Y->d = lda;
251: Y->xtype = CHOLMOD_SCALAR_TYPE;
252: Y->dtype = CHOLMOD_DOUBLE;
253: PetscFunctionReturn(PETSC_SUCCESS);
254: }
256: PetscErrorCode MatDenseUnWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y)
257: {
258: PetscFunctionBegin;
259: switch (rw) {
260: case GET_ARRAY_READ:
261: PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&Y->x));
262: break;
263: case GET_ARRAY_WRITE:
264: /* we don't have MatDenseRestoreArrayWrite */
265: PetscCall(MatDenseRestoreArray(X, (PetscScalar **)&Y->x));
266: break;
267: default:
268: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
269: break;
270: }
271: PetscFunctionReturn(PETSC_SUCCESS);
272: }
274: PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F)
275: {
276: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
278: PetscFunctionBegin;
279: if (chol->spqrfact) PetscCallExternal(!SuiteSparseQR_C_free, &chol->spqrfact, chol->common);
280: if (chol->factor) PetscCallExternal(!cholmod_X_free_factor, &chol->factor, chol->common);
281: if (chol->common->itype == CHOLMOD_INT) PetscCallExternal(!cholmod_finish, chol->common);
282: else PetscCallExternal(!cholmod_l_finish, chol->common);
283: PetscCall(PetscFree(chol->common));
284: PetscCall(PetscFree(chol->matrix));
285: PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", NULL));
286: PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorSymbolic_C", NULL));
287: PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", NULL));
288: PetscCall(PetscFree(F->data));
289: PetscFunctionReturn(PETSC_SUCCESS);
290: }
292: static PetscErrorCode MatSolve_CHOLMOD(Mat, Vec, Vec);
293: static PetscErrorCode MatMatSolve_CHOLMOD(Mat, Mat, Mat);
295: /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/
297: static PetscErrorCode MatView_Info_CHOLMOD(Mat F, PetscViewer viewer)
298: {
299: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
300: const cholmod_common *c = chol->common;
302: PetscFunctionBegin;
303: if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(PETSC_SUCCESS);
304: PetscCall(PetscViewerASCIIPrintf(viewer, "CHOLMOD run parameters:\n"));
305: PetscCall(PetscViewerASCIIPushTab(viewer));
306: PetscCall(PetscViewerASCIIPrintf(viewer, "Pack factors after symbolic factorization: %s\n", chol->pack ? "TRUE" : "FALSE"));
307: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.dbound %g (Smallest absolute value of diagonal entries of D)\n", c->dbound));
308: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow0 %g\n", c->grow0));
309: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow1 %g\n", c->grow1));
310: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow2 %u\n", (unsigned)c->grow2));
311: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.maxrank %u\n", (unsigned)c->maxrank));
312: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal_switch %g\n", c->supernodal_switch));
313: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal %d\n", c->supernodal));
314: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_asis %d\n", c->final_asis));
315: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_super %d\n", c->final_super));
316: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_ll %d\n", c->final_ll));
317: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_pack %d\n", c->final_pack));
318: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_monotonic %d\n", c->final_monotonic));
319: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_resymbol %d\n", c->final_resymbol));
320: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.zrelax [%g,%g,%g]\n", c->zrelax[0], c->zrelax[1], c->zrelax[2]));
321: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrelax [%u,%u,%u]\n", (unsigned)c->nrelax[0], (unsigned)c->nrelax[1], (unsigned)c->nrelax[2]));
322: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.prefer_upper %d\n", c->prefer_upper));
323: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.print %d\n", c->print));
324: for (PetscInt i = 0; i < c->nmethods; i++) {
325: PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering method %" PetscInt_FMT "%s:\n", i, i == c->selected ? " [SELECTED]" : ""));
326: PetscCall(PetscViewerASCIIPrintf(viewer, " lnz %g, fl %g, prune_dense %g, prune_dense2 %g\n", c->method[i].lnz, c->method[i].fl, c->method[i].prune_dense, c->method[i].prune_dense2));
327: }
328: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.postorder %d\n", c->postorder));
329: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n", c->default_nesdis));
330: /* Statistics */
331: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.fl %g (flop count from most recent analysis)\n", c->fl));
332: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.lnz %g (fundamental nz in L)\n", c->lnz));
333: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.anz %g\n", c->anz));
334: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.modfl %g (flop count from most recent update)\n", c->modfl));
335: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.malloc_count %g (number of live objects)\n", (double)c->malloc_count));
336: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_usage %g (peak memory usage in bytes)\n", (double)c->memory_usage));
337: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_inuse %g (current memory usage in bytes)\n", (double)c->memory_inuse));
338: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_col %g (number of column reallocations)\n", c->nrealloc_col));
339: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n", c->nrealloc_factor));
340: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n", c->ndbounds_hit));
341: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n", c->rowfacfl));
342: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n", c->aatfl));
343: #if defined(PETSC_USE_SUITESPARSE_GPU)
344: PetscCall(PetscViewerASCIIPrintf(viewer, "Common.useGPU %d\n", c->useGPU));
345: #endif
346: PetscCall(PetscViewerASCIIPopTab(viewer));
347: PetscFunctionReturn(PETSC_SUCCESS);
348: }
350: PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F, PetscViewer viewer)
351: {
352: PetscBool isascii;
353: PetscViewerFormat format;
355: PetscFunctionBegin;
356: PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
357: if (isascii) {
358: PetscCall(PetscViewerGetFormat(viewer, &format));
359: if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_CHOLMOD(F, viewer));
360: }
361: PetscFunctionReturn(PETSC_SUCCESS);
362: }
364: static PetscErrorCode MatSolve_CHOLMOD(Mat F, Vec B, Vec X)
365: {
366: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
367: cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL;
369: PetscFunctionBegin;
370: if (!F->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
371: static_F = F;
372: PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB));
373: PetscCall(VecWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
374: X_handle = &cholX;
375: PetscCallExternal(!cholmod_X_solve2, CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common);
376: PetscCallExternal(!cholmod_X_free_dense, &Y_handle, chol->common);
377: PetscCallExternal(!cholmod_X_free_dense, &E_handle, chol->common);
378: PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
379: PetscCall(VecUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
380: PetscCall(PetscLogFlops(4.0 * chol->common->lnz));
381: PetscFunctionReturn(PETSC_SUCCESS);
382: }
384: static PetscErrorCode MatMatSolve_CHOLMOD(Mat F, Mat B, Mat X)
385: {
386: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
387: cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL;
389: PetscFunctionBegin;
390: static_F = F;
391: PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB));
392: PetscCall(MatDenseWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
393: X_handle = &cholX;
394: PetscCallExternal(!cholmod_X_solve2, CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common);
395: PetscCallExternal(!cholmod_X_free_dense, &Y_handle, chol->common);
396: PetscCallExternal(!cholmod_X_free_dense, &E_handle, chol->common);
397: PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
398: PetscCall(MatDenseUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
399: PetscCall(PetscLogFlops(4.0 * B->cmap->n * chol->common->lnz));
400: PetscFunctionReturn(PETSC_SUCCESS);
401: }
403: static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F, Mat A, const MatFactorInfo *info)
404: {
405: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
406: cholmod_sparse cholA;
407: PetscBool aijalloc, valloc;
408: int err;
410: PetscFunctionBegin;
411: if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
412: PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc));
413: static_F = F;
414: err = !cholmod_X_factorize(&cholA, chol->factor, chol->common);
415: PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD factorization failed with status %d", chol->common->status);
416: PetscCheck(chol->common->status != CHOLMOD_NOT_POSDEF, PetscObjectComm((PetscObject)F), PETSC_ERR_MAT_CH_ZRPVT, "CHOLMOD detected that the matrix is not positive definite, failure at column %u", (unsigned)chol->factor->minor);
418: PetscCall(PetscLogFlops(chol->common->fl));
419: if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
420: if (valloc) PetscCall(PetscFree(cholA.x));
421: #if defined(PETSC_USE_SUITESPARSE_GPU)
422: PetscCall(PetscLogGpuTimeAdd(chol->common->CHOLMOD_GPU_GEMM_TIME + chol->common->CHOLMOD_GPU_SYRK_TIME + chol->common->CHOLMOD_GPU_TRSM_TIME + chol->common->CHOLMOD_GPU_POTRF_TIME));
423: #endif
425: F->ops->solve = MatSolve_CHOLMOD;
426: F->ops->solvetranspose = MatSolve_CHOLMOD;
427: F->ops->matsolve = MatMatSolve_CHOLMOD;
428: F->ops->matsolvetranspose = MatMatSolve_CHOLMOD;
429: PetscFunctionReturn(PETSC_SUCCESS);
430: }
432: PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F, Mat A, IS perm, const MatFactorInfo *info)
433: {
434: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
435: int err;
436: cholmod_sparse cholA;
437: PetscBool aijalloc, valloc;
438: PetscInt *fset = 0;
439: size_t fsize = 0;
441: PetscFunctionBegin;
442: F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD;
443: if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
444: /* Set options to F */
445: PetscCall(CholmodSetOptions(F));
447: PetscCall((*chol->Wrap)(A, PETSC_FALSE, &cholA, &aijalloc, &valloc));
448: static_F = F;
449: if (chol->factor) {
450: err = !cholmod_X_resymbol(&cholA, fset, fsize, (int)chol->pack, chol->factor, chol->common);
451: PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed with status %d", chol->common->status);
452: } else if (perm) {
453: const PetscInt *ip;
454: PetscCall(ISGetIndices(perm, &ip));
455: chol->factor = cholmod_X_analyze_p(&cholA, (PetscInt *)ip, fset, fsize, chol->common);
456: PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using PETSc ordering with status %d", chol->common->status);
457: PetscCall(ISRestoreIndices(perm, &ip));
458: } else {
459: chol->factor = cholmod_X_analyze(&cholA, chol->common);
460: PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status);
461: }
463: if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
464: if (valloc) PetscCall(PetscFree(cholA.x));
465: PetscFunctionReturn(PETSC_SUCCESS);
466: }
468: static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A, MatSolverType *type)
469: {
470: PetscFunctionBegin;
471: *type = MATSOLVERCHOLMOD;
472: PetscFunctionReturn(PETSC_SUCCESS);
473: }
475: PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F, MatInfoType flag, MatInfo *info)
476: {
477: Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
479: PetscFunctionBegin;
480: info->block_size = 1.0;
481: info->nz_allocated = chol->common->lnz;
482: info->nz_used = chol->common->lnz;
483: info->nz_unneeded = 0.0;
484: info->assemblies = 0.0;
485: info->mallocs = 0.0;
486: info->memory = chol->common->memory_inuse;
487: info->fill_ratio_given = 0;
488: info->fill_ratio_needed = 0;
489: info->factor_mallocs = chol->common->malloc_count;
490: PetscFunctionReturn(PETSC_SUCCESS);
491: }
493: /*MC
494: MATSOLVERCHOLMOD
496: A matrix type providing direct solvers (Cholesky) for sequential matrices
497: via the external package CHOLMOD.
499: Use `./configure --download-suitesparse` to install PETSc to use CHOLMOD
501: Use `-pc_type cholesky` `-pc_factor_mat_solver_type cholmod` to use this direct solver
503: Consult CHOLMOD documentation for more information about the common parameters
504: which correspond to the options database keys below.
506: Options Database Keys:
507: + -mat_cholmod_dbound 0 - Minimum absolute value of diagonal entries of D (None)
508: . -mat_cholmod_grow0 1.2 - Global growth ratio when factors are modified (None)
509: . -mat_cholmod_grow1 1.2 - Column growth ratio when factors are modified (None)
510: . -mat_cholmod_grow2 5 - Affine column growth constant when factors are modified (None)
511: . -mat_cholmod_maxrank 8 - Max rank of update, larger values are faster but use more memory [2,4,8] (None)
512: . -mat_cholmod_factor AUTO - (choose one of) `SIMPLICIAL`, `AUTO`, `SUPERNODAL`
513: . -mat_cholmod_supernodal_switch 40 - flop/nnz_L threshold for switching to supernodal factorization (None)
514: . -mat_cholmod_final_asis TRUE - Leave factors "as is" (None)
515: . -mat_cholmod_final_pack TRUE - Pack the columns when finished (use FALSE if the factors will be updated later) (None)
516: . -mat_cholmod_zrelax 0.8 - 3 real supernodal relaxed amalgamation parameters (None)
517: . -mat_cholmod_nrelax 4 - 3 size_t supernodal relaxed amalgamation parameters (None)
518: . -mat_cholmod_prefer_upper TRUE - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None)
519: . -mat_cholmod_print 3 - Verbosity level (None)
520: - -mat_ordering_type internal - Use the ordering provided by Cholmod
522: Level: beginner
524: Note:
525: CHOLMOD is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html>
527: .seealso: [](ch_matrices), `Mat`, `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType`
528: M*/
530: PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A, MatFactorType ftype, Mat *F)
531: {
532: Mat B;
533: Mat_CHOLMOD *chol;
534: PetscInt m = A->rmap->n, n = A->cmap->n, bs;
536: PetscFunctionBegin;
537: PetscCall(MatGetBlockSize(A, &bs));
538: *F = NULL;
539: if (bs != 1) {
540: PetscCall(PetscInfo(A, "CHOLMOD only supports block size=1.\n"));
541: PetscFunctionReturn(PETSC_SUCCESS);
542: }
543: if (PetscDefined(USE_COMPLEX) && A->hermitian != PETSC_BOOL3_TRUE) {
544: PetscCall(PetscInfo(A, "Only for Hermitian matrices.\n"));
545: PetscFunctionReturn(PETSC_SUCCESS);
546: }
547: /* Create the factorization matrix F */
548: PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
549: PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
550: PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name));
551: PetscCall(MatSetUp(B));
552: PetscCall(PetscNew(&chol));
554: chol->Wrap = MatWrapCholmod_seqsbaij;
555: B->data = chol;
557: B->ops->getinfo = MatGetInfo_CHOLMOD;
558: B->ops->view = MatView_CHOLMOD;
559: B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
560: B->ops->destroy = MatDestroy_CHOLMOD;
561: PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqsbaij_cholmod));
562: B->factortype = MAT_FACTOR_CHOLESKY;
563: B->assembled = PETSC_TRUE;
564: B->preallocated = PETSC_TRUE;
566: PetscCall(CholmodStart(B));
568: PetscCall(PetscFree(B->solvertype));
569: PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype));
570: B->canuseordering = PETSC_TRUE;
571: PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
572: *F = B;
573: PetscFunctionReturn(PETSC_SUCCESS);
574: }