1 2 /* 3 Provides an interface to the CHOLMOD sparse solver available through SuiteSparse version 4.2.1 4 5 When built with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the 6 integer type in UMFPACK, otherwise it will use int. This means 7 all integers in this file as simply declared as PetscInt. Also it means 8 that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit only] 9 10 */ 11 12 #include <../src/mat/impls/sbaij/seq/sbaij.h> 13 #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 14 15 /* 16 This is a terrible hack, but it allows the error handler to retain a context. 17 Note that this hack really cannot be made both reentrant and concurrent. 18 */ 19 static Mat static_F; 20 21 static void CholmodErrorHandler(int status, const char *file, int line, const char *message) { 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 } 32 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) 39 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) 46 47 #define CHOLMOD_OPTION_SIZE_T(name, help) \ 48 do { \ 49 PetscReal tmp = (PetscInt)c->name; \ 50 PetscCall(PetscOptionsReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 51 PetscCheck(tmp >= 0, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "value must be positive"); \ 52 c->name = (size_t)tmp; \ 53 } while (0) 54 55 #define CHOLMOD_OPTION_BOOL(name, help) \ 56 do { \ 57 PetscBool tmp = (PetscBool) !!c->name; \ 58 PetscCall(PetscOptionsBool("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 59 c->name = (int)tmp; \ 60 } while (0) 61 62 static PetscErrorCode CholmodSetOptions(Mat F) { 63 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 64 cholmod_common *c = chol->common; 65 PetscBool flg; 66 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"); 70 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 77 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; 82 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(0); 121 } 122 123 PetscErrorCode CholmodStart(Mat F) { 124 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 125 cholmod_common *c; 126 127 PetscFunctionBegin; 128 if (chol->common) PetscFunctionReturn(0); 129 PetscCall(PetscMalloc1(1, &chol->common)); 130 PetscCall(!cholmod_X_start(chol->common)); 131 132 c = chol->common; 133 c->error_handler = CholmodErrorHandler; 134 PetscFunctionReturn(0); 135 } 136 137 static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) { 138 Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ *)A->data; 139 PetscBool vallocin = PETSC_FALSE; 140 141 PetscFunctionBegin; 142 PetscCall(PetscMemzero(C, sizeof(*C))); 143 /* 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 */ 144 C->nrow = (size_t)A->cmap->n; 145 C->ncol = (size_t)A->rmap->n; 146 C->nzmax = (size_t)sbaij->maxnz; 147 C->p = sbaij->i; 148 C->i = sbaij->j; 149 if (values) { 150 #if defined(PETSC_USE_COMPLEX) 151 /* we need to pass CHOLMOD the conjugate matrix */ 152 PetscScalar *v; 153 PetscInt i; 154 155 PetscCall(PetscMalloc1(sbaij->maxnz, &v)); 156 for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]); 157 C->x = v; 158 vallocin = PETSC_TRUE; 159 #else 160 C->x = sbaij->a; 161 #endif 162 } 163 C->stype = -1; 164 C->itype = CHOLMOD_INT_TYPE; 165 C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 166 C->dtype = CHOLMOD_DOUBLE; 167 C->sorted = 1; 168 C->packed = 1; 169 *aijalloc = PETSC_FALSE; 170 *valloc = vallocin; 171 PetscFunctionReturn(0); 172 } 173 174 #define GET_ARRAY_READ 0 175 #define GET_ARRAY_WRITE 1 176 177 PetscErrorCode VecWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) { 178 PetscScalar *x; 179 PetscInt n; 180 181 PetscFunctionBegin; 182 PetscCall(PetscMemzero(Y, sizeof(*Y))); 183 switch (rw) { 184 case GET_ARRAY_READ: PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); break; 185 case GET_ARRAY_WRITE: PetscCall(VecGetArrayWrite(X, &x)); break; 186 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break; 187 } 188 PetscCall(VecGetSize(X, &n)); 189 190 Y->x = x; 191 Y->nrow = n; 192 Y->ncol = 1; 193 Y->nzmax = n; 194 Y->d = n; 195 Y->xtype = CHOLMOD_SCALAR_TYPE; 196 Y->dtype = CHOLMOD_DOUBLE; 197 PetscFunctionReturn(0); 198 } 199 200 PetscErrorCode VecUnWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) { 201 PetscFunctionBegin; 202 switch (rw) { 203 case GET_ARRAY_READ: PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&Y->x)); break; 204 case GET_ARRAY_WRITE: PetscCall(VecRestoreArrayWrite(X, (PetscScalar **)&Y->x)); break; 205 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break; 206 } 207 PetscFunctionReturn(0); 208 } 209 210 PetscErrorCode MatDenseWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) { 211 PetscScalar *x; 212 PetscInt m, n, lda; 213 214 PetscFunctionBegin; 215 PetscCall(PetscMemzero(Y, sizeof(*Y))); 216 switch (rw) { 217 case GET_ARRAY_READ: PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&x)); break; 218 case GET_ARRAY_WRITE: PetscCall(MatDenseGetArrayWrite(X, &x)); break; 219 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break; 220 } 221 PetscCall(MatDenseGetLDA(X, &lda)); 222 PetscCall(MatGetLocalSize(X, &m, &n)); 223 224 Y->x = x; 225 Y->nrow = m; 226 Y->ncol = n; 227 Y->nzmax = lda * n; 228 Y->d = lda; 229 Y->xtype = CHOLMOD_SCALAR_TYPE; 230 Y->dtype = CHOLMOD_DOUBLE; 231 PetscFunctionReturn(0); 232 } 233 234 PetscErrorCode MatDenseUnWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) { 235 PetscFunctionBegin; 236 switch (rw) { 237 case GET_ARRAY_READ: PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&Y->x)); break; 238 case GET_ARRAY_WRITE: 239 /* we don't have MatDenseRestoreArrayWrite */ 240 PetscCall(MatDenseRestoreArray(X, (PetscScalar **)&Y->x)); 241 break; 242 default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break; 243 } 244 PetscFunctionReturn(0); 245 } 246 247 PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) { 248 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 249 250 PetscFunctionBegin; 251 if (chol->spqrfact) PetscCall(!SuiteSparseQR_C_free(&chol->spqrfact, chol->common)); 252 if (chol->factor) PetscCall(!cholmod_X_free_factor(&chol->factor, chol->common)); 253 if (chol->common->itype == CHOLMOD_INT) { 254 PetscCall(!cholmod_finish(chol->common)); 255 } else { 256 PetscCall(!cholmod_l_finish(chol->common)); 257 } 258 PetscCall(PetscFree(chol->common)); 259 PetscCall(PetscFree(chol->matrix)); 260 PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", NULL)); 261 PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorSymbolic_C", NULL)); 262 PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", NULL)); 263 PetscCall(PetscFree(F->data)); 264 PetscFunctionReturn(0); 265 } 266 267 static PetscErrorCode MatSolve_CHOLMOD(Mat, Vec, Vec); 268 static PetscErrorCode MatMatSolve_CHOLMOD(Mat, Mat, Mat); 269 270 /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/ 271 272 static PetscErrorCode MatView_Info_CHOLMOD(Mat F, PetscViewer viewer) { 273 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 274 const cholmod_common *c = chol->common; 275 PetscInt i; 276 277 PetscFunctionBegin; 278 if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(0); 279 PetscCall(PetscViewerASCIIPrintf(viewer, "CHOLMOD run parameters:\n")); 280 PetscCall(PetscViewerASCIIPushTab(viewer)); 281 PetscCall(PetscViewerASCIIPrintf(viewer, "Pack factors after symbolic factorization: %s\n", chol->pack ? "TRUE" : "FALSE")); 282 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.dbound %g (Smallest absolute value of diagonal entries of D)\n", c->dbound)); 283 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow0 %g\n", c->grow0)); 284 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow1 %g\n", c->grow1)); 285 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow2 %u\n", (unsigned)c->grow2)); 286 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.maxrank %u\n", (unsigned)c->maxrank)); 287 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal_switch %g\n", c->supernodal_switch)); 288 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal %d\n", c->supernodal)); 289 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_asis %d\n", c->final_asis)); 290 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_super %d\n", c->final_super)); 291 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_ll %d\n", c->final_ll)); 292 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_pack %d\n", c->final_pack)); 293 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_monotonic %d\n", c->final_monotonic)); 294 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_resymbol %d\n", c->final_resymbol)); 295 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.zrelax [%g,%g,%g]\n", c->zrelax[0], c->zrelax[1], c->zrelax[2])); 296 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrelax [%u,%u,%u]\n", (unsigned)c->nrelax[0], (unsigned)c->nrelax[1], (unsigned)c->nrelax[2])); 297 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.prefer_upper %d\n", c->prefer_upper)); 298 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.print %d\n", c->print)); 299 for (i = 0; i < c->nmethods; i++) { 300 PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering method %" PetscInt_FMT "%s:\n", i, i == c->selected ? " [SELECTED]" : "")); 301 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)); 302 } 303 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.postorder %d\n", c->postorder)); 304 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n", c->default_nesdis)); 305 /* Statistics */ 306 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.fl %g (flop count from most recent analysis)\n", c->fl)); 307 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.lnz %g (fundamental nz in L)\n", c->lnz)); 308 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.anz %g\n", c->anz)); 309 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.modfl %g (flop count from most recent update)\n", c->modfl)); 310 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.malloc_count %g (number of live objects)\n", (double)c->malloc_count)); 311 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_usage %g (peak memory usage in bytes)\n", (double)c->memory_usage)); 312 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_inuse %g (current memory usage in bytes)\n", (double)c->memory_inuse)); 313 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_col %g (number of column reallocations)\n", c->nrealloc_col)); 314 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n", c->nrealloc_factor)); 315 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n", c->ndbounds_hit)); 316 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n", c->rowfacfl)); 317 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n", c->aatfl)); 318 #if defined(PETSC_USE_SUITESPARSE_GPU) 319 PetscCall(PetscViewerASCIIPrintf(viewer, "Common.useGPU %d\n", c->useGPU)); 320 #endif 321 PetscCall(PetscViewerASCIIPopTab(viewer)); 322 PetscFunctionReturn(0); 323 } 324 325 PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F, PetscViewer viewer) { 326 PetscBool iascii; 327 PetscViewerFormat format; 328 329 PetscFunctionBegin; 330 PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 331 if (iascii) { 332 PetscCall(PetscViewerGetFormat(viewer, &format)); 333 if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_CHOLMOD(F, viewer)); 334 } 335 PetscFunctionReturn(0); 336 } 337 338 static PetscErrorCode MatSolve_CHOLMOD(Mat F, Vec B, Vec X) { 339 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 340 cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL; 341 342 PetscFunctionBegin; 343 static_F = F; 344 PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB)); 345 PetscCall(VecWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 346 X_handle = &cholX; 347 PetscCall(!cholmod_X_solve2(CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common)); 348 PetscCall(!cholmod_X_free_dense(&Y_handle, chol->common)); 349 PetscCall(!cholmod_X_free_dense(&E_handle, chol->common)); 350 PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 351 PetscCall(VecUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 352 PetscCall(PetscLogFlops(4.0 * chol->common->lnz)); 353 PetscFunctionReturn(0); 354 } 355 356 static PetscErrorCode MatMatSolve_CHOLMOD(Mat F, Mat B, Mat X) { 357 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 358 cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL; 359 360 PetscFunctionBegin; 361 static_F = F; 362 PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB)); 363 PetscCall(MatDenseWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 364 X_handle = &cholX; 365 PetscCall(!cholmod_X_solve2(CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common)); 366 PetscCall(!cholmod_X_free_dense(&Y_handle, chol->common)); 367 PetscCall(!cholmod_X_free_dense(&E_handle, chol->common)); 368 PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 369 PetscCall(MatDenseUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 370 PetscCall(PetscLogFlops(4.0 * B->cmap->n * chol->common->lnz)); 371 PetscFunctionReturn(0); 372 } 373 374 static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F, Mat A, const MatFactorInfo *info) { 375 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 376 cholmod_sparse cholA; 377 PetscBool aijalloc, valloc; 378 int err; 379 380 PetscFunctionBegin; 381 PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc)); 382 static_F = F; 383 err = !cholmod_X_factorize(&cholA, chol->factor, chol->common); 384 PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD factorization failed with status %d", chol->common->status); 385 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); 386 387 PetscCall(PetscLogFlops(chol->common->fl)); 388 if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 389 if (valloc) PetscCall(PetscFree(cholA.x)); 390 #if defined(PETSC_USE_SUITESPARSE_GPU) 391 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)); 392 #endif 393 394 F->ops->solve = MatSolve_CHOLMOD; 395 F->ops->solvetranspose = MatSolve_CHOLMOD; 396 F->ops->matsolve = MatMatSolve_CHOLMOD; 397 F->ops->matsolvetranspose = MatMatSolve_CHOLMOD; 398 PetscFunctionReturn(0); 399 } 400 401 PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F, Mat A, IS perm, const MatFactorInfo *info) { 402 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 403 int err; 404 cholmod_sparse cholA; 405 PetscBool aijalloc, valloc; 406 PetscInt *fset = 0; 407 size_t fsize = 0; 408 409 PetscFunctionBegin; 410 /* Set options to F */ 411 PetscCall(CholmodSetOptions(F)); 412 413 PetscCall((*chol->Wrap)(A, PETSC_FALSE, &cholA, &aijalloc, &valloc)); 414 static_F = F; 415 if (chol->factor) { 416 err = !cholmod_X_resymbol(&cholA, fset, fsize, (int)chol->pack, chol->factor, chol->common); 417 PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed with status %d", chol->common->status); 418 } else if (perm) { 419 const PetscInt *ip; 420 PetscCall(ISGetIndices(perm, &ip)); 421 chol->factor = cholmod_X_analyze_p(&cholA, (PetscInt *)ip, fset, fsize, chol->common); 422 PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using PETSc ordering with status %d", chol->common->status); 423 PetscCall(ISRestoreIndices(perm, &ip)); 424 } else { 425 chol->factor = cholmod_X_analyze(&cholA, chol->common); 426 PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status); 427 } 428 429 if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 430 if (valloc) PetscCall(PetscFree(cholA.x)); 431 432 F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD; 433 PetscFunctionReturn(0); 434 } 435 436 static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A, MatSolverType *type) { 437 PetscFunctionBegin; 438 *type = MATSOLVERCHOLMOD; 439 PetscFunctionReturn(0); 440 } 441 442 PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F, MatInfoType flag, MatInfo *info) { 443 Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 444 445 PetscFunctionBegin; 446 info->block_size = 1.0; 447 info->nz_allocated = chol->common->lnz; 448 info->nz_used = chol->common->lnz; 449 info->nz_unneeded = 0.0; 450 info->assemblies = 0.0; 451 info->mallocs = 0.0; 452 info->memory = chol->common->memory_inuse; 453 info->fill_ratio_given = 0; 454 info->fill_ratio_needed = 0; 455 info->factor_mallocs = chol->common->malloc_count; 456 PetscFunctionReturn(0); 457 } 458 459 /*MC 460 MATSOLVERCHOLMOD 461 462 A matrix type providing direct solvers (Cholesky) for sequential matrices 463 via the external package CHOLMOD. 464 465 Use ./configure --download-suitesparse to install PETSc to use CHOLMOD 466 467 Use -pc_type cholesky -pc_factor_mat_solver_type cholmod to use this direct solver 468 469 Consult CHOLMOD documentation for more information about the Common parameters 470 which correspond to the options database keys below. 471 472 Options Database Keys: 473 + -mat_cholmod_dbound <0> - Minimum absolute value of diagonal entries of D (None) 474 . -mat_cholmod_grow0 <1.2> - Global growth ratio when factors are modified (None) 475 . -mat_cholmod_grow1 <1.2> - Column growth ratio when factors are modified (None) 476 . -mat_cholmod_grow2 <5> - Affine column growth constant when factors are modified (None) 477 . -mat_cholmod_maxrank <8> - Max rank of update, larger values are faster but use more memory [2,4,8] (None) 478 . -mat_cholmod_factor <AUTO> - (choose one of) SIMPLICIAL AUTO SUPERNODAL 479 . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None) 480 . -mat_cholmod_final_asis <TRUE> - Leave factors "as is" (None) 481 . -mat_cholmod_final_pack <TRUE> - Pack the columns when finished (use FALSE if the factors will be updated later) (None) 482 . -mat_cholmod_zrelax <0.8> - 3 real supernodal relaxed amalgamation parameters (None) 483 . -mat_cholmod_nrelax <4> - 3 size_t supernodal relaxed amalgamation parameters (None) 484 . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None) 485 . -mat_cholmod_print <3> - Verbosity level (None) 486 - -mat_ordering_type internal - Use the ordering provided by Cholmod 487 488 Level: beginner 489 490 Note: CHOLMOD is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 491 492 .seealso: `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType` 493 M*/ 494 495 PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A, MatFactorType ftype, Mat *F) { 496 Mat B; 497 Mat_CHOLMOD *chol; 498 PetscInt m = A->rmap->n, n = A->cmap->n, bs; 499 500 PetscFunctionBegin; 501 PetscCall(MatGetBlockSize(A, &bs)); 502 PetscCheck(bs == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "CHOLMOD only supports block size=1, given %" PetscInt_FMT, bs); 503 #if defined(PETSC_USE_COMPLEX) 504 PetscCheck(A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only for Hermitian matrices"); 505 #endif 506 /* Create the factorization matrix F */ 507 PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 508 PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 509 PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name)); 510 PetscCall(MatSetUp(B)); 511 PetscCall(PetscNew(&chol)); 512 513 chol->Wrap = MatWrapCholmod_seqsbaij; 514 B->data = chol; 515 516 B->ops->getinfo = MatGetInfo_CHOLMOD; 517 B->ops->view = MatView_CHOLMOD; 518 B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 519 B->ops->destroy = MatDestroy_CHOLMOD; 520 PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqsbaij_cholmod)); 521 B->factortype = MAT_FACTOR_CHOLESKY; 522 B->assembled = PETSC_TRUE; 523 B->preallocated = PETSC_TRUE; 524 525 PetscCall(CholmodStart(B)); 526 527 PetscCall(PetscFree(B->solvertype)); 528 PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype)); 529 B->canuseordering = PETSC_TRUE; 530 PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 531 *F = B; 532 PetscFunctionReturn(0); 533 } 534