1641875f9SMatthew G Knepley /* 2d515b9b4SShri Abhyankar Provides an interface to the CHOLMOD sparse solver available through SuiteSparse version 4.2.1 3641875f9SMatthew G Knepley 48999bf53SRichard Mills When built with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the 5641875f9SMatthew G Knepley integer type in UMFPACK, otherwise it will use int. This means 6641875f9SMatthew G Knepley all integers in this file as simply declared as PetscInt. Also it means 77de69702SBarry Smith that one cannot use 64BIT_INDICES on 32-bit pointer systems [as Suitesparse_long is 32-bit only] 8641875f9SMatthew G Knepley 9641875f9SMatthew G Knepley */ 10641875f9SMatthew G Knepley 11c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 12c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 13641875f9SMatthew G Knepley 14641875f9SMatthew G Knepley /* 15641875f9SMatthew G Knepley This is a terrible hack, but it allows the error handler to retain a context. 16641875f9SMatthew G Knepley Note that this hack really cannot be made both reentrant and concurrent. 17641875f9SMatthew G Knepley */ 18641875f9SMatthew G Knepley static Mat static_F; 19641875f9SMatthew G Knepley 20d71ae5a4SJacob Faibussowitsch static void CholmodErrorHandler(int status, const char *file, int line, const char *message) 21d71ae5a4SJacob Faibussowitsch { 22641875f9SMatthew G Knepley PetscFunctionBegin; 23641875f9SMatthew G Knepley if (status > CHOLMOD_OK) { 249566063dSJacob Faibussowitsch PetscCallVoid(PetscInfo(static_F, "CHOLMOD warning %d at %s:%d: %s\n", status, file, line, message)); 25641875f9SMatthew G Knepley } else if (status == CHOLMOD_OK) { /* Documentation says this can happen, but why? */ 269566063dSJacob Faibussowitsch PetscCallVoid(PetscInfo(static_F, "CHOLMOD OK at %s:%d: %s\n", file, line, message)); 27641875f9SMatthew G Knepley } else { 289566063dSJacob Faibussowitsch PetscCallVoid(PetscErrorPrintf("CHOLMOD error %d at %s:%d: %s\n", status, file, line, message)); 29641875f9SMatthew G Knepley } 30641875f9SMatthew G Knepley PetscFunctionReturnVoid(); 31641875f9SMatthew G Knepley } 32641875f9SMatthew G Knepley 339371c9d4SSatish Balay #define CHOLMOD_OPTION_DOUBLE(name, help) \ 349371c9d4SSatish Balay do { \ 35641875f9SMatthew G Knepley PetscReal tmp = (PetscReal)c->name; \ 369566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 37641875f9SMatthew G Knepley c->name = (double)tmp; \ 38641875f9SMatthew G Knepley } while (0) 3926fbe8dcSKarl Rupp 409371c9d4SSatish Balay #define CHOLMOD_OPTION_INT(name, help) \ 419371c9d4SSatish Balay do { \ 42641875f9SMatthew G Knepley PetscInt tmp = (PetscInt)c->name; \ 439566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 44641875f9SMatthew G Knepley c->name = (int)tmp; \ 45641875f9SMatthew G Knepley } while (0) 4626fbe8dcSKarl Rupp 479371c9d4SSatish Balay #define CHOLMOD_OPTION_SIZE_T(name, help) \ 489371c9d4SSatish Balay do { \ 4954b3d318SStefano Zampini PetscReal tmp = (PetscInt)c->name; \ 509566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 5108401ef6SPierre Jolivet PetscCheck(tmp >= 0, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "value must be positive"); \ 52641875f9SMatthew G Knepley c->name = (size_t)tmp; \ 53641875f9SMatthew G Knepley } while (0) 5426fbe8dcSKarl Rupp 559371c9d4SSatish Balay #define CHOLMOD_OPTION_BOOL(name, help) \ 569371c9d4SSatish Balay do { \ 57ace3abfcSBarry Smith PetscBool tmp = (PetscBool) !!c->name; \ 589566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 59641875f9SMatthew G Knepley c->name = (int)tmp; \ 60641875f9SMatthew G Knepley } while (0) 61641875f9SMatthew G Knepley 62d71ae5a4SJacob Faibussowitsch static PetscErrorCode CholmodSetOptions(Mat F) 63d71ae5a4SJacob Faibussowitsch { 6426cc229bSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 6526cc229bSBarry Smith cholmod_common *c = chol->common; 6626cc229bSBarry Smith PetscBool flg; 6726cc229bSBarry Smith 6826cc229bSBarry Smith PetscFunctionBegin; 69d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "CHOLMOD Options", "Mat"); 7054b3d318SStefano Zampini CHOLMOD_OPTION_INT(nmethods, "Number of different ordering methods to try"); 7126fbe8dcSKarl Rupp 72b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 73b9eaa5e8SBarry Smith c->useGPU = 1; 74b9eaa5e8SBarry Smith CHOLMOD_OPTION_INT(useGPU, "Use GPU for BLAS 1, otherwise 0"); 7554b3d318SStefano Zampini CHOLMOD_OPTION_SIZE_T(maxGpuMemBytes, "Maximum memory to allocate on the GPU"); 7654b3d318SStefano Zampini CHOLMOD_OPTION_DOUBLE(maxGpuMemFraction, "Fraction of available GPU memory to allocate"); 77b9eaa5e8SBarry Smith #endif 78b9eaa5e8SBarry Smith 7954b3d318SStefano Zampini /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */ 8054b3d318SStefano Zampini chol->pack = (PetscBool)c->final_pack; 819566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_pack", "Pack factors after factorization [disable for frequent repeat factorization]", "None", chol->pack, &chol->pack, NULL)); 82641875f9SMatthew G Knepley c->final_pack = (int)chol->pack; 83641875f9SMatthew G Knepley 84641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(dbound, "Minimum absolute value of diagonal entries of D"); 85641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow0, "Global growth ratio when factors are modified"); 86641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow1, "Column growth ratio when factors are modified"); 87641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(grow2, "Affine column growth constant when factors are modified"); 88641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(maxrank, "Max rank of update, larger values are faster but use more memory [2,4,8]"); 89641875f9SMatthew G Knepley { 90641875f9SMatthew G Knepley static const char *const list[] = {"SIMPLICIAL", "AUTO", "SUPERNODAL", "MatCholmodFactorType", "MAT_CHOLMOD_FACTOR_", 0}; 919566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_cholmod_factor", "Factorization method", "None", list, (PetscEnum)c->supernodal, (PetscEnum *)&c->supernodal, NULL)); 92641875f9SMatthew G Knepley } 93641875f9SMatthew G Knepley if (c->supernodal) CHOLMOD_OPTION_DOUBLE(supernodal_switch, "flop/nnz_L threshold for switching to supernodal factorization"); 94b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_asis, "Leave factors \"as is\""); 95b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_pack, "Pack the columns when finished (use FALSE if the factors will be updated later)"); 96641875f9SMatthew G Knepley if (!c->final_asis) { 97b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_super, "Leave supernodal factors instead of converting to simplicial"); 98b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_ll, "Turn LDL' factorization into LL'"); 99b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_monotonic, "Ensure columns are monotonic when done"); 100b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_resymbol, "Remove numerically zero values resulting from relaxed supernodal amalgamation"); 101641875f9SMatthew G Knepley } 102641875f9SMatthew G Knepley { 103641875f9SMatthew G Knepley PetscReal tmp[] = {(PetscReal)c->zrelax[0], (PetscReal)c->zrelax[1], (PetscReal)c->zrelax[2]}; 104641875f9SMatthew G Knepley PetscInt n = 3; 1059566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-mat_cholmod_zrelax", "3 real supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg)); 106aed4548fSBarry Smith PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_zrelax"); 1079371c9d4SSatish Balay if (flg) 1089371c9d4SSatish Balay while (n--) c->zrelax[n] = (double)tmp[n]; 109641875f9SMatthew G Knepley } 110641875f9SMatthew G Knepley { 111641875f9SMatthew G Knepley PetscInt n, tmp[] = {(PetscInt)c->nrelax[0], (PetscInt)c->nrelax[1], (PetscInt)c->nrelax[2]}; 1129566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-mat_cholmod_nrelax", "3 size_t supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg)); 113aed4548fSBarry Smith PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_nrelax"); 1149371c9d4SSatish Balay if (flg) 1159371c9d4SSatish Balay while (n--) c->nrelax[n] = (size_t)tmp[n]; 116641875f9SMatthew G Knepley } 117b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(prefer_upper, "Work with upper triangular form [faster when using fill-reducing ordering, slower in natural ordering]"); 118b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(default_nesdis, "Use NESDIS instead of METIS for nested dissection"); 119641875f9SMatthew G Knepley CHOLMOD_OPTION_INT(print, "Verbosity level"); 120d0609cedSBarry Smith PetscOptionsEnd(); 1213ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 122641875f9SMatthew G Knepley } 123641875f9SMatthew G Knepley 124d71ae5a4SJacob Faibussowitsch PetscErrorCode CholmodStart(Mat F) 125d71ae5a4SJacob Faibussowitsch { 12626cc229bSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 12726cc229bSBarry Smith cholmod_common *c; 12826cc229bSBarry Smith 12926cc229bSBarry Smith PetscFunctionBegin; 1303ba16761SJacob Faibussowitsch if (chol->common) PetscFunctionReturn(PETSC_SUCCESS); 13126cc229bSBarry Smith PetscCall(PetscMalloc1(1, &chol->common)); 1323ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_start, chol->common); 13326cc229bSBarry Smith 13426cc229bSBarry Smith c = chol->common; 13526cc229bSBarry Smith c->error_handler = CholmodErrorHandler; 1363ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13726cc229bSBarry Smith } 13826cc229bSBarry Smith 139d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) 140d71ae5a4SJacob Faibussowitsch { 141641875f9SMatthew G Knepley Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ *)A->data; 142218db3c1SStefano Zampini PetscBool vallocin = PETSC_FALSE; 143641875f9SMatthew G Knepley 144641875f9SMatthew G Knepley PetscFunctionBegin; 1459566063dSJacob Faibussowitsch PetscCall(PetscMemzero(C, sizeof(*C))); 146641875f9SMatthew G Knepley /* 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 */ 147641875f9SMatthew G Knepley C->nrow = (size_t)A->cmap->n; 148641875f9SMatthew G Knepley C->ncol = (size_t)A->rmap->n; 149641875f9SMatthew G Knepley C->nzmax = (size_t)sbaij->maxnz; 150641875f9SMatthew G Knepley C->p = sbaij->i; 151641875f9SMatthew G Knepley C->i = sbaij->j; 152218db3c1SStefano Zampini if (values) { 153218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 154218db3c1SStefano Zampini /* we need to pass CHOLMOD the conjugate matrix */ 155218db3c1SStefano Zampini PetscScalar *v; 156218db3c1SStefano Zampini PetscInt i; 157218db3c1SStefano Zampini 1589566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sbaij->maxnz, &v)); 159218db3c1SStefano Zampini for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]); 160218db3c1SStefano Zampini C->x = v; 161218db3c1SStefano Zampini vallocin = PETSC_TRUE; 162218db3c1SStefano Zampini #else 163641875f9SMatthew G Knepley C->x = sbaij->a; 164218db3c1SStefano Zampini #endif 165218db3c1SStefano Zampini } 166641875f9SMatthew G Knepley C->stype = -1; 167641875f9SMatthew G Knepley C->itype = CHOLMOD_INT_TYPE; 168218db3c1SStefano Zampini C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 169641875f9SMatthew G Knepley C->dtype = CHOLMOD_DOUBLE; 170641875f9SMatthew G Knepley C->sorted = 1; 171641875f9SMatthew G Knepley C->packed = 1; 172641875f9SMatthew G Knepley *aijalloc = PETSC_FALSE; 173218db3c1SStefano Zampini *valloc = vallocin; 1743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 175641875f9SMatthew G Knepley } 176641875f9SMatthew G Knepley 177218db3c1SStefano Zampini #define GET_ARRAY_READ 0 178218db3c1SStefano Zampini #define GET_ARRAY_WRITE 1 179218db3c1SStefano Zampini 180d71ae5a4SJacob Faibussowitsch PetscErrorCode VecWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) 181d71ae5a4SJacob Faibussowitsch { 182218db3c1SStefano Zampini PetscScalar *x; 183641875f9SMatthew G Knepley PetscInt n; 184641875f9SMatthew G Knepley 185641875f9SMatthew G Knepley PetscFunctionBegin; 1869566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y, sizeof(*Y))); 187218db3c1SStefano Zampini switch (rw) { 188d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 189d71ae5a4SJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 190d71ae5a4SJacob Faibussowitsch break; 191d71ae5a4SJacob Faibussowitsch case GET_ARRAY_WRITE: 192d71ae5a4SJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &x)); 193d71ae5a4SJacob Faibussowitsch break; 194d71ae5a4SJacob Faibussowitsch default: 195d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 196d71ae5a4SJacob Faibussowitsch break; 197218db3c1SStefano Zampini } 1989566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &n)); 19926fbe8dcSKarl Rupp 200218db3c1SStefano Zampini Y->x = x; 201641875f9SMatthew G Knepley Y->nrow = n; 202641875f9SMatthew G Knepley Y->ncol = 1; 203641875f9SMatthew G Knepley Y->nzmax = n; 204641875f9SMatthew G Knepley Y->d = n; 205641875f9SMatthew G Knepley Y->xtype = CHOLMOD_SCALAR_TYPE; 206641875f9SMatthew G Knepley Y->dtype = CHOLMOD_DOUBLE; 2073ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 208641875f9SMatthew G Knepley } 209641875f9SMatthew G Knepley 210d71ae5a4SJacob Faibussowitsch PetscErrorCode VecUnWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) 211d71ae5a4SJacob Faibussowitsch { 212d9ca1df4SBarry Smith PetscFunctionBegin; 213218db3c1SStefano Zampini switch (rw) { 214d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 215d71ae5a4SJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&Y->x)); 216d71ae5a4SJacob Faibussowitsch break; 217d71ae5a4SJacob Faibussowitsch case GET_ARRAY_WRITE: 218d71ae5a4SJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, (PetscScalar **)&Y->x)); 219d71ae5a4SJacob Faibussowitsch break; 220d71ae5a4SJacob Faibussowitsch default: 221d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 222d71ae5a4SJacob Faibussowitsch break; 223218db3c1SStefano Zampini } 2243ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 225218db3c1SStefano Zampini } 226218db3c1SStefano Zampini 227d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDenseWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) 228d71ae5a4SJacob Faibussowitsch { 229218db3c1SStefano Zampini PetscScalar *x; 230218db3c1SStefano Zampini PetscInt m, n, lda; 231218db3c1SStefano Zampini 232218db3c1SStefano Zampini PetscFunctionBegin; 2339566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y, sizeof(*Y))); 234218db3c1SStefano Zampini switch (rw) { 235d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 236d71ae5a4SJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&x)); 237d71ae5a4SJacob Faibussowitsch break; 238d71ae5a4SJacob Faibussowitsch case GET_ARRAY_WRITE: 239d71ae5a4SJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(X, &x)); 240d71ae5a4SJacob Faibussowitsch break; 241d71ae5a4SJacob Faibussowitsch default: 242d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 243d71ae5a4SJacob Faibussowitsch break; 244218db3c1SStefano Zampini } 2459566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X, &lda)); 2469566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m, &n)); 247218db3c1SStefano Zampini 248218db3c1SStefano Zampini Y->x = x; 249218db3c1SStefano Zampini Y->nrow = m; 250218db3c1SStefano Zampini Y->ncol = n; 251218db3c1SStefano Zampini Y->nzmax = lda * n; 252218db3c1SStefano Zampini Y->d = lda; 253218db3c1SStefano Zampini Y->xtype = CHOLMOD_SCALAR_TYPE; 254218db3c1SStefano Zampini Y->dtype = CHOLMOD_DOUBLE; 2553ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 256218db3c1SStefano Zampini } 257218db3c1SStefano Zampini 258d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDenseUnWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) 259d71ae5a4SJacob Faibussowitsch { 260218db3c1SStefano Zampini PetscFunctionBegin; 261218db3c1SStefano Zampini switch (rw) { 262d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 263d71ae5a4SJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&Y->x)); 264d71ae5a4SJacob Faibussowitsch break; 265218db3c1SStefano Zampini case GET_ARRAY_WRITE: 266218db3c1SStefano Zampini /* we don't have MatDenseRestoreArrayWrite */ 2679566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(X, (PetscScalar **)&Y->x)); 268218db3c1SStefano Zampini break; 269d71ae5a4SJacob Faibussowitsch default: 270d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 271d71ae5a4SJacob Faibussowitsch break; 272218db3c1SStefano Zampini } 2733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 274d9ca1df4SBarry Smith } 275d9ca1df4SBarry Smith 276d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) 277d71ae5a4SJacob Faibussowitsch { 2786b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 279641875f9SMatthew G Knepley 280641875f9SMatthew G Knepley PetscFunctionBegin; 2813ba16761SJacob Faibussowitsch if (chol->spqrfact) PetscCallExternal(!SuiteSparseQR_C_free, &chol->spqrfact, chol->common); 2823ba16761SJacob Faibussowitsch if (chol->factor) PetscCallExternal(!cholmod_X_free_factor, &chol->factor, chol->common); 283a2fc1e05SToby Isaac if (chol->common->itype == CHOLMOD_INT) { 2843ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_finish, chol->common); 285a2fc1e05SToby Isaac } else { 2863ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_finish, chol->common); 287a2fc1e05SToby Isaac } 2889566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->common)); 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->matrix)); 2909566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", NULL)); 2912e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorSymbolic_C", NULL)); 2922e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", NULL)); 2939566063dSJacob Faibussowitsch PetscCall(PetscFree(F->data)); 2943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 295641875f9SMatthew G Knepley } 296641875f9SMatthew G Knepley 297641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat, Vec, Vec); 298218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat, Mat, Mat); 299641875f9SMatthew G Knepley 300fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/ 301641875f9SMatthew G Knepley 302d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_CHOLMOD(Mat F, PetscViewer viewer) 303d71ae5a4SJacob Faibussowitsch { 3046b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 305641875f9SMatthew G Knepley const cholmod_common *c = chol->common; 306641875f9SMatthew G Knepley PetscInt i; 307641875f9SMatthew G Knepley 308641875f9SMatthew G Knepley PetscFunctionBegin; 3093ba16761SJacob Faibussowitsch if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(PETSC_SUCCESS); 3109566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CHOLMOD run parameters:\n")); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Pack factors after symbolic factorization: %s\n", chol->pack ? "TRUE" : "FALSE")); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.dbound %g (Smallest absolute value of diagonal entries of D)\n", c->dbound)); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow0 %g\n", c->grow0)); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow1 %g\n", c->grow1)); 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow2 %u\n", (unsigned)c->grow2)); 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.maxrank %u\n", (unsigned)c->maxrank)); 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal_switch %g\n", c->supernodal_switch)); 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal %d\n", c->supernodal)); 3209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_asis %d\n", c->final_asis)); 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_super %d\n", c->final_super)); 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_ll %d\n", c->final_ll)); 3239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_pack %d\n", c->final_pack)); 3249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_monotonic %d\n", c->final_monotonic)); 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_resymbol %d\n", c->final_resymbol)); 3269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.zrelax [%g,%g,%g]\n", c->zrelax[0], c->zrelax[1], c->zrelax[2])); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrelax [%u,%u,%u]\n", (unsigned)c->nrelax[0], (unsigned)c->nrelax[1], (unsigned)c->nrelax[2])); 3289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.prefer_upper %d\n", c->prefer_upper)); 3299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.print %d\n", c->print)); 330641875f9SMatthew G Knepley for (i = 0; i < c->nmethods; i++) { 3319566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering method %" PetscInt_FMT "%s:\n", i, i == c->selected ? " [SELECTED]" : "")); 3329371c9d4SSatish Balay 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)); 333641875f9SMatthew G Knepley } 3349566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.postorder %d\n", c->postorder)); 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n", c->default_nesdis)); 336641875f9SMatthew G Knepley /* Statistics */ 3379566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.fl %g (flop count from most recent analysis)\n", c->fl)); 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.lnz %g (fundamental nz in L)\n", c->lnz)); 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.anz %g\n", c->anz)); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.modfl %g (flop count from most recent update)\n", c->modfl)); 3419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.malloc_count %g (number of live objects)\n", (double)c->malloc_count)); 3429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_usage %g (peak memory usage in bytes)\n", (double)c->memory_usage)); 3439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_inuse %g (current memory usage in bytes)\n", (double)c->memory_inuse)); 3449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_col %g (number of column reallocations)\n", c->nrealloc_col)); 3459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n", c->nrealloc_factor)); 3469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n", c->ndbounds_hit)); 3479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n", c->rowfacfl)); 3489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n", c->aatfl)); 349b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 3509566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.useGPU %d\n", c->useGPU)); 351b9eaa5e8SBarry Smith #endif 3529566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3533ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 354641875f9SMatthew G Knepley } 355641875f9SMatthew G Knepley 356d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F, PetscViewer viewer) 357d71ae5a4SJacob Faibussowitsch { 358ace3abfcSBarry Smith PetscBool iascii; 359641875f9SMatthew G Knepley PetscViewerFormat format; 360641875f9SMatthew G Knepley 361641875f9SMatthew G Knepley PetscFunctionBegin; 3629566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 363641875f9SMatthew G Knepley if (iascii) { 3649566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 36548a46eb9SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_CHOLMOD(F, viewer)); 366641875f9SMatthew G Knepley } 3673ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 368641875f9SMatthew G Knepley } 369641875f9SMatthew G Knepley 370d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_CHOLMOD(Mat F, Vec B, Vec X) 371d71ae5a4SJacob Faibussowitsch { 3726b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 373218db3c1SStefano Zampini cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL; 374641875f9SMatthew G Knepley 375641875f9SMatthew G Knepley PetscFunctionBegin; 376641875f9SMatthew G Knepley static_F = F; 3779566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB)); 3789566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 379218db3c1SStefano Zampini X_handle = &cholX; 3803ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_solve2, CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common); 3813ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &Y_handle, chol->common); 3823ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &E_handle, chol->common); 3839566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 3849566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 3859566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * chol->common->lnz)); 3863ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 387218db3c1SStefano Zampini } 388218db3c1SStefano Zampini 389d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_CHOLMOD(Mat F, Mat B, Mat X) 390d71ae5a4SJacob Faibussowitsch { 391218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 392218db3c1SStefano Zampini cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL; 393218db3c1SStefano Zampini 394218db3c1SStefano Zampini PetscFunctionBegin; 395218db3c1SStefano Zampini static_F = F; 3969566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB)); 3979566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 398218db3c1SStefano Zampini X_handle = &cholX; 3993ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_solve2, CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common); 4003ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &Y_handle, chol->common); 4013ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &E_handle, chol->common); 4029566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 4039566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 4049566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * B->cmap->n * chol->common->lnz)); 4053ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 406641875f9SMatthew G Knepley } 407641875f9SMatthew G Knepley 408d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F, Mat A, const MatFactorInfo *info) 409d71ae5a4SJacob Faibussowitsch { 4106b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 411641875f9SMatthew G Knepley cholmod_sparse cholA; 412218db3c1SStefano Zampini PetscBool aijalloc, valloc; 413d0609cedSBarry Smith int err; 414641875f9SMatthew G Knepley 415641875f9SMatthew G Knepley PetscFunctionBegin; 4169566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc)); 417641875f9SMatthew G Knepley static_F = F; 418d0609cedSBarry Smith err = !cholmod_X_factorize(&cholA, chol->factor, chol->common); 419d0609cedSBarry Smith PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD factorization failed with status %d", chol->common->status); 42008401ef6SPierre Jolivet 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); 421641875f9SMatthew G Knepley 4229566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(chol->common->fl)); 4239566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 4249566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 425218db3c1SStefano Zampini #if defined(PETSC_USE_SUITESPARSE_GPU) 4269566063dSJacob Faibussowitsch 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)); 427218db3c1SStefano Zampini #endif 428641875f9SMatthew G Knepley 429641875f9SMatthew G Knepley F->ops->solve = MatSolve_CHOLMOD; 430641875f9SMatthew G Knepley F->ops->solvetranspose = MatSolve_CHOLMOD; 431218db3c1SStefano Zampini F->ops->matsolve = MatMatSolve_CHOLMOD; 432218db3c1SStefano Zampini F->ops->matsolvetranspose = MatMatSolve_CHOLMOD; 4333ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 434641875f9SMatthew G Knepley } 435641875f9SMatthew G Knepley 436d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F, Mat A, IS perm, const MatFactorInfo *info) 437d71ae5a4SJacob Faibussowitsch { 4386b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 439d0609cedSBarry Smith int err; 440641875f9SMatthew G Knepley cholmod_sparse cholA; 441218db3c1SStefano Zampini PetscBool aijalloc, valloc; 442641875f9SMatthew G Knepley PetscInt *fset = 0; 443641875f9SMatthew G Knepley size_t fsize = 0; 444641875f9SMatthew G Knepley 445641875f9SMatthew G Knepley PetscFunctionBegin; 44626cc229bSBarry Smith /* Set options to F */ 44726cc229bSBarry Smith PetscCall(CholmodSetOptions(F)); 44826cc229bSBarry Smith 4499566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A, PETSC_FALSE, &cholA, &aijalloc, &valloc)); 450641875f9SMatthew G Knepley static_F = F; 451641875f9SMatthew G Knepley if (chol->factor) { 452d0609cedSBarry Smith err = !cholmod_X_resymbol(&cholA, fset, fsize, (int)chol->pack, chol->factor, chol->common); 453d0609cedSBarry Smith PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed with status %d", chol->common->status); 454641875f9SMatthew G Knepley } else if (perm) { 455641875f9SMatthew G Knepley const PetscInt *ip; 4569566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &ip)); 457641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze_p(&cholA, (PetscInt *)ip, fset, fsize, chol->common); 45828b400f6SJacob Faibussowitsch PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using PETSc ordering with status %d", chol->common->status); 4599566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &ip)); 460641875f9SMatthew G Knepley } else { 461641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze(&cholA, chol->common); 46228b400f6SJacob Faibussowitsch PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status); 463641875f9SMatthew G Knepley } 464641875f9SMatthew G Knepley 4659566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 4669566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 467641875f9SMatthew G Knepley 468641875f9SMatthew G Knepley F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD; 4693ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 470641875f9SMatthew G Knepley } 471641875f9SMatthew G Knepley 472d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A, MatSolverType *type) 473d71ae5a4SJacob Faibussowitsch { 474641875f9SMatthew G Knepley PetscFunctionBegin; 475641875f9SMatthew G Knepley *type = MATSOLVERCHOLMOD; 4763ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 477641875f9SMatthew G Knepley } 478641875f9SMatthew G Knepley 479d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F, MatInfoType flag, MatInfo *info) 480d71ae5a4SJacob Faibussowitsch { 481218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 482218db3c1SStefano Zampini 483218db3c1SStefano Zampini PetscFunctionBegin; 484218db3c1SStefano Zampini info->block_size = 1.0; 485218db3c1SStefano Zampini info->nz_allocated = chol->common->lnz; 486218db3c1SStefano Zampini info->nz_used = chol->common->lnz; 487218db3c1SStefano Zampini info->nz_unneeded = 0.0; 488218db3c1SStefano Zampini info->assemblies = 0.0; 489218db3c1SStefano Zampini info->mallocs = 0.0; 490218db3c1SStefano Zampini info->memory = chol->common->memory_inuse; 491218db3c1SStefano Zampini info->fill_ratio_given = 0; 492218db3c1SStefano Zampini info->fill_ratio_needed = 0; 493218db3c1SStefano Zampini info->factor_mallocs = chol->common->malloc_count; 4943ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 495218db3c1SStefano Zampini } 496218db3c1SStefano Zampini 497641875f9SMatthew G Knepley /*MC 498ee300463SSatish Balay MATSOLVERCHOLMOD 499ee300463SSatish Balay 500ee300463SSatish Balay A matrix type providing direct solvers (Cholesky) for sequential matrices 501641875f9SMatthew G Knepley via the external package CHOLMOD. 502641875f9SMatthew G Knepley 5032ef1f0ffSBarry Smith Use `./configure --download-suitesparse` to install PETSc to use CHOLMOD 504c2b89b5dSBarry Smith 5052ef1f0ffSBarry Smith Use `-pc_type cholesky` `-pc_factor_mat_solver_type cholmod` to use this direct solver 506641875f9SMatthew G Knepley 5072ef1f0ffSBarry Smith Consult CHOLMOD documentation for more information about the common parameters 508641875f9SMatthew G Knepley which correspond to the options database keys below. 509641875f9SMatthew G Knepley 510641875f9SMatthew G Knepley Options Database Keys: 511e08999f5SMatthew G Knepley + -mat_cholmod_dbound <0> - Minimum absolute value of diagonal entries of D (None) 512e08999f5SMatthew G Knepley . -mat_cholmod_grow0 <1.2> - Global growth ratio when factors are modified (None) 513e08999f5SMatthew G Knepley . -mat_cholmod_grow1 <1.2> - Column growth ratio when factors are modified (None) 514e08999f5SMatthew G Knepley . -mat_cholmod_grow2 <5> - Affine column growth constant when factors are modified (None) 515e08999f5SMatthew G Knepley . -mat_cholmod_maxrank <8> - Max rank of update, larger values are faster but use more memory [2,4,8] (None) 5162ef1f0ffSBarry Smith . -mat_cholmod_factor <AUTO> - (choose one of) `SIMPLICIAL`, `AUTO`, `SUPERNODAL` 517e08999f5SMatthew G Knepley . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None) 518e08999f5SMatthew G Knepley . -mat_cholmod_final_asis <TRUE> - Leave factors "as is" (None) 519e08999f5SMatthew G Knepley . -mat_cholmod_final_pack <TRUE> - Pack the columns when finished (use FALSE if the factors will be updated later) (None) 520e08999f5SMatthew G Knepley . -mat_cholmod_zrelax <0.8> - 3 real supernodal relaxed amalgamation parameters (None) 521e08999f5SMatthew G Knepley . -mat_cholmod_nrelax <4> - 3 size_t supernodal relaxed amalgamation parameters (None) 522e08999f5SMatthew G Knepley . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None) 5232c7c0729SBarry Smith . -mat_cholmod_print <3> - Verbosity level (None) 524ee300463SSatish Balay - -mat_ordering_type internal - Use the ordering provided by Cholmod 525641875f9SMatthew G Knepley 526641875f9SMatthew G Knepley Level: beginner 527641875f9SMatthew G Knepley 5282ef1f0ffSBarry Smith Note: 529*1d27aa22SBarry Smith CHOLMOD is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html> 530a364b7d2SBarry Smith 5311cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType` 532641875f9SMatthew G Knepley M*/ 533b2573a8aSBarry Smith 534d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A, MatFactorType ftype, Mat *F) 535d71ae5a4SJacob Faibussowitsch { 536641875f9SMatthew G Knepley Mat B; 537641875f9SMatthew G Knepley Mat_CHOLMOD *chol; 538641875f9SMatthew G Knepley PetscInt m = A->rmap->n, n = A->cmap->n, bs; 539641875f9SMatthew G Knepley 540641875f9SMatthew G Knepley PetscFunctionBegin; 5419566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 54203e5aca4SStefano Zampini *F = NULL; 54303e5aca4SStefano Zampini if (bs != 1) { 54403e5aca4SStefano Zampini PetscCall(PetscInfo(A, "CHOLMOD only supports block size=1.\n")); 54503e5aca4SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 54603e5aca4SStefano Zampini } 547218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 54803e5aca4SStefano Zampini if (A->hermitian != PETSC_BOOL3_TRUE) { 54903e5aca4SStefano Zampini PetscCall(PetscInfo(A, "Only for Hermitian matrices.\n")); 55003e5aca4SStefano Zampini PetscFunctionReturn(PETSC_SUCCESS); 55103e5aca4SStefano Zampini } 552218db3c1SStefano Zampini #endif 553641875f9SMatthew G Knepley /* Create the factorization matrix F */ 5549566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 5559566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 5569566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name)); 5579566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 5584dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&chol)); 55926fbe8dcSKarl Rupp 560641875f9SMatthew G Knepley chol->Wrap = MatWrapCholmod_seqsbaij; 5616b8f6f9dSBarry Smith B->data = chol; 562641875f9SMatthew G Knepley 563218db3c1SStefano Zampini B->ops->getinfo = MatGetInfo_CHOLMOD; 564641875f9SMatthew G Knepley B->ops->view = MatView_CHOLMOD; 565641875f9SMatthew G Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 566641875f9SMatthew G Knepley B->ops->destroy = MatDestroy_CHOLMOD; 5679566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqsbaij_cholmod)); 568641875f9SMatthew G Knepley B->factortype = MAT_FACTOR_CHOLESKY; 569218db3c1SStefano Zampini B->assembled = PETSC_TRUE; 570641875f9SMatthew G Knepley B->preallocated = PETSC_TRUE; 571641875f9SMatthew G Knepley 5729566063dSJacob Faibussowitsch PetscCall(CholmodStart(B)); 57300c67f3bSHong Zhang 5749566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 5759566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype)); 576f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 5779566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 578641875f9SMatthew G Knepley *F = B; 5793ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 580641875f9SMatthew G Knepley } 581