1641875f9SMatthew G Knepley 2641875f9SMatthew G Knepley /* 3d515b9b4SShri Abhyankar Provides an interface to the CHOLMOD sparse solver available through SuiteSparse version 4.2.1 4641875f9SMatthew G Knepley 58999bf53SRichard Mills When built with PETSC_USE_64BIT_INDICES this will use Suitesparse_long as the 6641875f9SMatthew G Knepley integer type in UMFPACK, otherwise it will use int. This means 7641875f9SMatthew G Knepley all integers in this file as simply declared as PetscInt. Also it means 8*7de69702SBarry Smith that one cannot use 64BIT_INDICES on 32-bit pointer systems [as Suitesparse_long is 32-bit only] 9641875f9SMatthew G Knepley 10641875f9SMatthew G Knepley */ 11641875f9SMatthew G Knepley 12c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/sbaij.h> 13c6db04a5SJed Brown #include <../src/mat/impls/sbaij/seq/cholmod/cholmodimpl.h> 14641875f9SMatthew G Knepley 15641875f9SMatthew G Knepley /* 16641875f9SMatthew G Knepley This is a terrible hack, but it allows the error handler to retain a context. 17641875f9SMatthew G Knepley Note that this hack really cannot be made both reentrant and concurrent. 18641875f9SMatthew G Knepley */ 19641875f9SMatthew G Knepley static Mat static_F; 20641875f9SMatthew G Knepley 21d71ae5a4SJacob Faibussowitsch static void CholmodErrorHandler(int status, const char *file, int line, const char *message) 22d71ae5a4SJacob Faibussowitsch { 23641875f9SMatthew G Knepley PetscFunctionBegin; 24641875f9SMatthew G Knepley if (status > CHOLMOD_OK) { 259566063dSJacob Faibussowitsch PetscCallVoid(PetscInfo(static_F, "CHOLMOD warning %d at %s:%d: %s\n", status, file, line, message)); 26641875f9SMatthew G Knepley } else if (status == CHOLMOD_OK) { /* Documentation says this can happen, but why? */ 279566063dSJacob Faibussowitsch PetscCallVoid(PetscInfo(static_F, "CHOLMOD OK at %s:%d: %s\n", file, line, message)); 28641875f9SMatthew G Knepley } else { 299566063dSJacob Faibussowitsch PetscCallVoid(PetscErrorPrintf("CHOLMOD error %d at %s:%d: %s\n", status, file, line, message)); 30641875f9SMatthew G Knepley } 31641875f9SMatthew G Knepley PetscFunctionReturnVoid(); 32641875f9SMatthew G Knepley } 33641875f9SMatthew G Knepley 349371c9d4SSatish Balay #define CHOLMOD_OPTION_DOUBLE(name, help) \ 359371c9d4SSatish Balay do { \ 36641875f9SMatthew G Knepley PetscReal tmp = (PetscReal)c->name; \ 379566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 38641875f9SMatthew G Knepley c->name = (double)tmp; \ 39641875f9SMatthew G Knepley } while (0) 4026fbe8dcSKarl Rupp 419371c9d4SSatish Balay #define CHOLMOD_OPTION_INT(name, help) \ 429371c9d4SSatish Balay do { \ 43641875f9SMatthew G Knepley PetscInt tmp = (PetscInt)c->name; \ 449566063dSJacob Faibussowitsch PetscCall(PetscOptionsInt("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 45641875f9SMatthew G Knepley c->name = (int)tmp; \ 46641875f9SMatthew G Knepley } while (0) 4726fbe8dcSKarl Rupp 489371c9d4SSatish Balay #define CHOLMOD_OPTION_SIZE_T(name, help) \ 499371c9d4SSatish Balay do { \ 5054b3d318SStefano Zampini PetscReal tmp = (PetscInt)c->name; \ 519566063dSJacob Faibussowitsch PetscCall(PetscOptionsReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 5208401ef6SPierre Jolivet PetscCheck(tmp >= 0, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "value must be positive"); \ 53641875f9SMatthew G Knepley c->name = (size_t)tmp; \ 54641875f9SMatthew G Knepley } while (0) 5526fbe8dcSKarl Rupp 569371c9d4SSatish Balay #define CHOLMOD_OPTION_BOOL(name, help) \ 579371c9d4SSatish Balay do { \ 58ace3abfcSBarry Smith PetscBool tmp = (PetscBool) !!c->name; \ 599566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \ 60641875f9SMatthew G Knepley c->name = (int)tmp; \ 61641875f9SMatthew G Knepley } while (0) 62641875f9SMatthew G Knepley 63d71ae5a4SJacob Faibussowitsch static PetscErrorCode CholmodSetOptions(Mat F) 64d71ae5a4SJacob Faibussowitsch { 6526cc229bSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 6626cc229bSBarry Smith cholmod_common *c = chol->common; 6726cc229bSBarry Smith PetscBool flg; 6826cc229bSBarry Smith 6926cc229bSBarry Smith PetscFunctionBegin; 70d0609cedSBarry Smith PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "CHOLMOD Options", "Mat"); 7154b3d318SStefano Zampini CHOLMOD_OPTION_INT(nmethods, "Number of different ordering methods to try"); 7226fbe8dcSKarl Rupp 73b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 74b9eaa5e8SBarry Smith c->useGPU = 1; 75b9eaa5e8SBarry Smith CHOLMOD_OPTION_INT(useGPU, "Use GPU for BLAS 1, otherwise 0"); 7654b3d318SStefano Zampini CHOLMOD_OPTION_SIZE_T(maxGpuMemBytes, "Maximum memory to allocate on the GPU"); 7754b3d318SStefano Zampini CHOLMOD_OPTION_DOUBLE(maxGpuMemFraction, "Fraction of available GPU memory to allocate"); 78b9eaa5e8SBarry Smith #endif 79b9eaa5e8SBarry Smith 8054b3d318SStefano Zampini /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */ 8154b3d318SStefano Zampini chol->pack = (PetscBool)c->final_pack; 829566063dSJacob Faibussowitsch PetscCall(PetscOptionsBool("-mat_cholmod_pack", "Pack factors after factorization [disable for frequent repeat factorization]", "None", chol->pack, &chol->pack, NULL)); 83641875f9SMatthew G Knepley c->final_pack = (int)chol->pack; 84641875f9SMatthew G Knepley 85641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(dbound, "Minimum absolute value of diagonal entries of D"); 86641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow0, "Global growth ratio when factors are modified"); 87641875f9SMatthew G Knepley CHOLMOD_OPTION_DOUBLE(grow1, "Column growth ratio when factors are modified"); 88641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(grow2, "Affine column growth constant when factors are modified"); 89641875f9SMatthew G Knepley CHOLMOD_OPTION_SIZE_T(maxrank, "Max rank of update, larger values are faster but use more memory [2,4,8]"); 90641875f9SMatthew G Knepley { 91641875f9SMatthew G Knepley static const char *const list[] = {"SIMPLICIAL", "AUTO", "SUPERNODAL", "MatCholmodFactorType", "MAT_CHOLMOD_FACTOR_", 0}; 929566063dSJacob Faibussowitsch PetscCall(PetscOptionsEnum("-mat_cholmod_factor", "Factorization method", "None", list, (PetscEnum)c->supernodal, (PetscEnum *)&c->supernodal, NULL)); 93641875f9SMatthew G Knepley } 94641875f9SMatthew G Knepley if (c->supernodal) CHOLMOD_OPTION_DOUBLE(supernodal_switch, "flop/nnz_L threshold for switching to supernodal factorization"); 95b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_asis, "Leave factors \"as is\""); 96b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_pack, "Pack the columns when finished (use FALSE if the factors will be updated later)"); 97641875f9SMatthew G Knepley if (!c->final_asis) { 98b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_super, "Leave supernodal factors instead of converting to simplicial"); 99b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_ll, "Turn LDL' factorization into LL'"); 100b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_monotonic, "Ensure columns are monotonic when done"); 101b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(final_resymbol, "Remove numerically zero values resulting from relaxed supernodal amalgamation"); 102641875f9SMatthew G Knepley } 103641875f9SMatthew G Knepley { 104641875f9SMatthew G Knepley PetscReal tmp[] = {(PetscReal)c->zrelax[0], (PetscReal)c->zrelax[1], (PetscReal)c->zrelax[2]}; 105641875f9SMatthew G Knepley PetscInt n = 3; 1069566063dSJacob Faibussowitsch PetscCall(PetscOptionsRealArray("-mat_cholmod_zrelax", "3 real supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg)); 107aed4548fSBarry Smith PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_zrelax"); 1089371c9d4SSatish Balay if (flg) 1099371c9d4SSatish Balay while (n--) c->zrelax[n] = (double)tmp[n]; 110641875f9SMatthew G Knepley } 111641875f9SMatthew G Knepley { 112641875f9SMatthew G Knepley PetscInt n, tmp[] = {(PetscInt)c->nrelax[0], (PetscInt)c->nrelax[1], (PetscInt)c->nrelax[2]}; 1139566063dSJacob Faibussowitsch PetscCall(PetscOptionsIntArray("-mat_cholmod_nrelax", "3 size_t supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg)); 114aed4548fSBarry Smith PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_nrelax"); 1159371c9d4SSatish Balay if (flg) 1169371c9d4SSatish Balay while (n--) c->nrelax[n] = (size_t)tmp[n]; 117641875f9SMatthew G Knepley } 118b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(prefer_upper, "Work with upper triangular form [faster when using fill-reducing ordering, slower in natural ordering]"); 119b9eaa5e8SBarry Smith CHOLMOD_OPTION_BOOL(default_nesdis, "Use NESDIS instead of METIS for nested dissection"); 120641875f9SMatthew G Knepley CHOLMOD_OPTION_INT(print, "Verbosity level"); 121d0609cedSBarry Smith PetscOptionsEnd(); 1223ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 123641875f9SMatthew G Knepley } 124641875f9SMatthew G Knepley 125d71ae5a4SJacob Faibussowitsch PetscErrorCode CholmodStart(Mat F) 126d71ae5a4SJacob Faibussowitsch { 12726cc229bSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 12826cc229bSBarry Smith cholmod_common *c; 12926cc229bSBarry Smith 13026cc229bSBarry Smith PetscFunctionBegin; 1313ba16761SJacob Faibussowitsch if (chol->common) PetscFunctionReturn(PETSC_SUCCESS); 13226cc229bSBarry Smith PetscCall(PetscMalloc1(1, &chol->common)); 1333ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_start, chol->common); 13426cc229bSBarry Smith 13526cc229bSBarry Smith c = chol->common; 13626cc229bSBarry Smith c->error_handler = CholmodErrorHandler; 1373ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 13826cc229bSBarry Smith } 13926cc229bSBarry Smith 140d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) 141d71ae5a4SJacob Faibussowitsch { 142641875f9SMatthew G Knepley Mat_SeqSBAIJ *sbaij = (Mat_SeqSBAIJ *)A->data; 143218db3c1SStefano Zampini PetscBool vallocin = PETSC_FALSE; 144641875f9SMatthew G Knepley 145641875f9SMatthew G Knepley PetscFunctionBegin; 1469566063dSJacob Faibussowitsch PetscCall(PetscMemzero(C, sizeof(*C))); 147641875f9SMatthew 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 */ 148641875f9SMatthew G Knepley C->nrow = (size_t)A->cmap->n; 149641875f9SMatthew G Knepley C->ncol = (size_t)A->rmap->n; 150641875f9SMatthew G Knepley C->nzmax = (size_t)sbaij->maxnz; 151641875f9SMatthew G Knepley C->p = sbaij->i; 152641875f9SMatthew G Knepley C->i = sbaij->j; 153218db3c1SStefano Zampini if (values) { 154218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 155218db3c1SStefano Zampini /* we need to pass CHOLMOD the conjugate matrix */ 156218db3c1SStefano Zampini PetscScalar *v; 157218db3c1SStefano Zampini PetscInt i; 158218db3c1SStefano Zampini 1599566063dSJacob Faibussowitsch PetscCall(PetscMalloc1(sbaij->maxnz, &v)); 160218db3c1SStefano Zampini for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]); 161218db3c1SStefano Zampini C->x = v; 162218db3c1SStefano Zampini vallocin = PETSC_TRUE; 163218db3c1SStefano Zampini #else 164641875f9SMatthew G Knepley C->x = sbaij->a; 165218db3c1SStefano Zampini #endif 166218db3c1SStefano Zampini } 167641875f9SMatthew G Knepley C->stype = -1; 168641875f9SMatthew G Knepley C->itype = CHOLMOD_INT_TYPE; 169218db3c1SStefano Zampini C->xtype = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN; 170641875f9SMatthew G Knepley C->dtype = CHOLMOD_DOUBLE; 171641875f9SMatthew G Knepley C->sorted = 1; 172641875f9SMatthew G Knepley C->packed = 1; 173641875f9SMatthew G Knepley *aijalloc = PETSC_FALSE; 174218db3c1SStefano Zampini *valloc = vallocin; 1753ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 176641875f9SMatthew G Knepley } 177641875f9SMatthew G Knepley 178218db3c1SStefano Zampini #define GET_ARRAY_READ 0 179218db3c1SStefano Zampini #define GET_ARRAY_WRITE 1 180218db3c1SStefano Zampini 181d71ae5a4SJacob Faibussowitsch PetscErrorCode VecWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) 182d71ae5a4SJacob Faibussowitsch { 183218db3c1SStefano Zampini PetscScalar *x; 184641875f9SMatthew G Knepley PetscInt n; 185641875f9SMatthew G Knepley 186641875f9SMatthew G Knepley PetscFunctionBegin; 1879566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y, sizeof(*Y))); 188218db3c1SStefano Zampini switch (rw) { 189d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 190d71ae5a4SJacob Faibussowitsch PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); 191d71ae5a4SJacob Faibussowitsch break; 192d71ae5a4SJacob Faibussowitsch case GET_ARRAY_WRITE: 193d71ae5a4SJacob Faibussowitsch PetscCall(VecGetArrayWrite(X, &x)); 194d71ae5a4SJacob Faibussowitsch break; 195d71ae5a4SJacob Faibussowitsch default: 196d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 197d71ae5a4SJacob Faibussowitsch break; 198218db3c1SStefano Zampini } 1999566063dSJacob Faibussowitsch PetscCall(VecGetSize(X, &n)); 20026fbe8dcSKarl Rupp 201218db3c1SStefano Zampini Y->x = x; 202641875f9SMatthew G Knepley Y->nrow = n; 203641875f9SMatthew G Knepley Y->ncol = 1; 204641875f9SMatthew G Knepley Y->nzmax = n; 205641875f9SMatthew G Knepley Y->d = n; 206641875f9SMatthew G Knepley Y->xtype = CHOLMOD_SCALAR_TYPE; 207641875f9SMatthew G Knepley Y->dtype = CHOLMOD_DOUBLE; 2083ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 209641875f9SMatthew G Knepley } 210641875f9SMatthew G Knepley 211d71ae5a4SJacob Faibussowitsch PetscErrorCode VecUnWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) 212d71ae5a4SJacob Faibussowitsch { 213d9ca1df4SBarry Smith PetscFunctionBegin; 214218db3c1SStefano Zampini switch (rw) { 215d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 216d71ae5a4SJacob Faibussowitsch PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&Y->x)); 217d71ae5a4SJacob Faibussowitsch break; 218d71ae5a4SJacob Faibussowitsch case GET_ARRAY_WRITE: 219d71ae5a4SJacob Faibussowitsch PetscCall(VecRestoreArrayWrite(X, (PetscScalar **)&Y->x)); 220d71ae5a4SJacob Faibussowitsch break; 221d71ae5a4SJacob Faibussowitsch default: 222d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 223d71ae5a4SJacob Faibussowitsch break; 224218db3c1SStefano Zampini } 2253ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 226218db3c1SStefano Zampini } 227218db3c1SStefano Zampini 228d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDenseWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) 229d71ae5a4SJacob Faibussowitsch { 230218db3c1SStefano Zampini PetscScalar *x; 231218db3c1SStefano Zampini PetscInt m, n, lda; 232218db3c1SStefano Zampini 233218db3c1SStefano Zampini PetscFunctionBegin; 2349566063dSJacob Faibussowitsch PetscCall(PetscMemzero(Y, sizeof(*Y))); 235218db3c1SStefano Zampini switch (rw) { 236d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 237d71ae5a4SJacob Faibussowitsch PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&x)); 238d71ae5a4SJacob Faibussowitsch break; 239d71ae5a4SJacob Faibussowitsch case GET_ARRAY_WRITE: 240d71ae5a4SJacob Faibussowitsch PetscCall(MatDenseGetArrayWrite(X, &x)); 241d71ae5a4SJacob Faibussowitsch break; 242d71ae5a4SJacob Faibussowitsch default: 243d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 244d71ae5a4SJacob Faibussowitsch break; 245218db3c1SStefano Zampini } 2469566063dSJacob Faibussowitsch PetscCall(MatDenseGetLDA(X, &lda)); 2479566063dSJacob Faibussowitsch PetscCall(MatGetLocalSize(X, &m, &n)); 248218db3c1SStefano Zampini 249218db3c1SStefano Zampini Y->x = x; 250218db3c1SStefano Zampini Y->nrow = m; 251218db3c1SStefano Zampini Y->ncol = n; 252218db3c1SStefano Zampini Y->nzmax = lda * n; 253218db3c1SStefano Zampini Y->d = lda; 254218db3c1SStefano Zampini Y->xtype = CHOLMOD_SCALAR_TYPE; 255218db3c1SStefano Zampini Y->dtype = CHOLMOD_DOUBLE; 2563ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 257218db3c1SStefano Zampini } 258218db3c1SStefano Zampini 259d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDenseUnWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) 260d71ae5a4SJacob Faibussowitsch { 261218db3c1SStefano Zampini PetscFunctionBegin; 262218db3c1SStefano Zampini switch (rw) { 263d71ae5a4SJacob Faibussowitsch case GET_ARRAY_READ: 264d71ae5a4SJacob Faibussowitsch PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&Y->x)); 265d71ae5a4SJacob Faibussowitsch break; 266218db3c1SStefano Zampini case GET_ARRAY_WRITE: 267218db3c1SStefano Zampini /* we don't have MatDenseRestoreArrayWrite */ 2689566063dSJacob Faibussowitsch PetscCall(MatDenseRestoreArray(X, (PetscScalar **)&Y->x)); 269218db3c1SStefano Zampini break; 270d71ae5a4SJacob Faibussowitsch default: 271d71ae5a4SJacob Faibussowitsch SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); 272d71ae5a4SJacob Faibussowitsch break; 273218db3c1SStefano Zampini } 2743ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 275d9ca1df4SBarry Smith } 276d9ca1df4SBarry Smith 277d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) 278d71ae5a4SJacob Faibussowitsch { 2796b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 280641875f9SMatthew G Knepley 281641875f9SMatthew G Knepley PetscFunctionBegin; 2823ba16761SJacob Faibussowitsch if (chol->spqrfact) PetscCallExternal(!SuiteSparseQR_C_free, &chol->spqrfact, chol->common); 2833ba16761SJacob Faibussowitsch if (chol->factor) PetscCallExternal(!cholmod_X_free_factor, &chol->factor, chol->common); 284a2fc1e05SToby Isaac if (chol->common->itype == CHOLMOD_INT) { 2853ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_finish, chol->common); 286a2fc1e05SToby Isaac } else { 2873ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_l_finish, chol->common); 288a2fc1e05SToby Isaac } 2899566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->common)); 2909566063dSJacob Faibussowitsch PetscCall(PetscFree(chol->matrix)); 2919566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", NULL)); 2922e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorSymbolic_C", NULL)); 2932e956fe4SStefano Zampini PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", NULL)); 2949566063dSJacob Faibussowitsch PetscCall(PetscFree(F->data)); 2953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 296641875f9SMatthew G Knepley } 297641875f9SMatthew G Knepley 298641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat, Vec, Vec); 299218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat, Mat, Mat); 300641875f9SMatthew G Knepley 301fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/ 302641875f9SMatthew G Knepley 303d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_CHOLMOD(Mat F, PetscViewer viewer) 304d71ae5a4SJacob Faibussowitsch { 3056b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 306641875f9SMatthew G Knepley const cholmod_common *c = chol->common; 307641875f9SMatthew G Knepley PetscInt i; 308641875f9SMatthew G Knepley 309641875f9SMatthew G Knepley PetscFunctionBegin; 3103ba16761SJacob Faibussowitsch if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(PETSC_SUCCESS); 3119566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "CHOLMOD run parameters:\n")); 3129566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPushTab(viewer)); 3139566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Pack factors after symbolic factorization: %s\n", chol->pack ? "TRUE" : "FALSE")); 3149566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.dbound %g (Smallest absolute value of diagonal entries of D)\n", c->dbound)); 3159566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow0 %g\n", c->grow0)); 3169566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow1 %g\n", c->grow1)); 3179566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow2 %u\n", (unsigned)c->grow2)); 3189566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.maxrank %u\n", (unsigned)c->maxrank)); 3199566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal_switch %g\n", c->supernodal_switch)); 3209566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal %d\n", c->supernodal)); 3219566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_asis %d\n", c->final_asis)); 3229566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_super %d\n", c->final_super)); 3239566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_ll %d\n", c->final_ll)); 3249566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_pack %d\n", c->final_pack)); 3259566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_monotonic %d\n", c->final_monotonic)); 3269566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_resymbol %d\n", c->final_resymbol)); 3279566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.zrelax [%g,%g,%g]\n", c->zrelax[0], c->zrelax[1], c->zrelax[2])); 3289566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrelax [%u,%u,%u]\n", (unsigned)c->nrelax[0], (unsigned)c->nrelax[1], (unsigned)c->nrelax[2])); 3299566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.prefer_upper %d\n", c->prefer_upper)); 3309566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.print %d\n", c->print)); 331641875f9SMatthew G Knepley for (i = 0; i < c->nmethods; i++) { 3329566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering method %" PetscInt_FMT "%s:\n", i, i == c->selected ? " [SELECTED]" : "")); 3339371c9d4SSatish 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)); 334641875f9SMatthew G Knepley } 3359566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.postorder %d\n", c->postorder)); 3369566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.default_nesdis %d (use NESDIS instead of METIS for nested dissection)\n", c->default_nesdis)); 337641875f9SMatthew G Knepley /* Statistics */ 3389566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.fl %g (flop count from most recent analysis)\n", c->fl)); 3399566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.lnz %g (fundamental nz in L)\n", c->lnz)); 3409566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.anz %g\n", c->anz)); 3419566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.modfl %g (flop count from most recent update)\n", c->modfl)); 3429566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.malloc_count %g (number of live objects)\n", (double)c->malloc_count)); 3439566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_usage %g (peak memory usage in bytes)\n", (double)c->memory_usage)); 3449566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_inuse %g (current memory usage in bytes)\n", (double)c->memory_inuse)); 3459566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_col %g (number of column reallocations)\n", c->nrealloc_col)); 3469566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_factor %g (number of factor reallocations due to column reallocations)\n", c->nrealloc_factor)); 3479566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.ndbounds_hit %g (number of times diagonal was modified by dbound)\n", c->ndbounds_hit)); 3489566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.rowfacfl %g (number of flops in last call to cholmod_rowfac)\n", c->rowfacfl)); 3499566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.aatfl %g (number of flops to compute A(:,f)*A(:,f)')\n", c->aatfl)); 350b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU) 3519566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPrintf(viewer, "Common.useGPU %d\n", c->useGPU)); 352b9eaa5e8SBarry Smith #endif 3539566063dSJacob Faibussowitsch PetscCall(PetscViewerASCIIPopTab(viewer)); 3543ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 355641875f9SMatthew G Knepley } 356641875f9SMatthew G Knepley 357d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F, PetscViewer viewer) 358d71ae5a4SJacob Faibussowitsch { 359ace3abfcSBarry Smith PetscBool iascii; 360641875f9SMatthew G Knepley PetscViewerFormat format; 361641875f9SMatthew G Knepley 362641875f9SMatthew G Knepley PetscFunctionBegin; 3639566063dSJacob Faibussowitsch PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii)); 364641875f9SMatthew G Knepley if (iascii) { 3659566063dSJacob Faibussowitsch PetscCall(PetscViewerGetFormat(viewer, &format)); 36648a46eb9SPierre Jolivet if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_CHOLMOD(F, viewer)); 367641875f9SMatthew G Knepley } 3683ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 369641875f9SMatthew G Knepley } 370641875f9SMatthew G Knepley 371d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_CHOLMOD(Mat F, Vec B, Vec X) 372d71ae5a4SJacob Faibussowitsch { 3736b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 374218db3c1SStefano Zampini cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL; 375641875f9SMatthew G Knepley 376641875f9SMatthew G Knepley PetscFunctionBegin; 377641875f9SMatthew G Knepley static_F = F; 3789566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB)); 3799566063dSJacob Faibussowitsch PetscCall(VecWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 380218db3c1SStefano Zampini X_handle = &cholX; 3813ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_solve2, CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common); 3823ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &Y_handle, chol->common); 3833ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &E_handle, chol->common); 3849566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 3859566063dSJacob Faibussowitsch PetscCall(VecUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 3869566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * chol->common->lnz)); 3873ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 388218db3c1SStefano Zampini } 389218db3c1SStefano Zampini 390d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatMatSolve_CHOLMOD(Mat F, Mat B, Mat X) 391d71ae5a4SJacob Faibussowitsch { 392218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 393218db3c1SStefano Zampini cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL; 394218db3c1SStefano Zampini 395218db3c1SStefano Zampini PetscFunctionBegin; 396218db3c1SStefano Zampini static_F = F; 3979566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB)); 3989566063dSJacob Faibussowitsch PetscCall(MatDenseWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 399218db3c1SStefano Zampini X_handle = &cholX; 4003ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_solve2, CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common); 4013ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &Y_handle, chol->common); 4023ba16761SJacob Faibussowitsch PetscCallExternal(!cholmod_X_free_dense, &E_handle, chol->common); 4039566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB)); 4049566063dSJacob Faibussowitsch PetscCall(MatDenseUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX)); 4059566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(4.0 * B->cmap->n * chol->common->lnz)); 4063ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 407641875f9SMatthew G Knepley } 408641875f9SMatthew G Knepley 409d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F, Mat A, const MatFactorInfo *info) 410d71ae5a4SJacob Faibussowitsch { 4116b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 412641875f9SMatthew G Knepley cholmod_sparse cholA; 413218db3c1SStefano Zampini PetscBool aijalloc, valloc; 414d0609cedSBarry Smith int err; 415641875f9SMatthew G Knepley 416641875f9SMatthew G Knepley PetscFunctionBegin; 4179566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc)); 418641875f9SMatthew G Knepley static_F = F; 419d0609cedSBarry Smith err = !cholmod_X_factorize(&cholA, chol->factor, chol->common); 420d0609cedSBarry Smith PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD factorization failed with status %d", chol->common->status); 42108401ef6SPierre 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); 422641875f9SMatthew G Knepley 4239566063dSJacob Faibussowitsch PetscCall(PetscLogFlops(chol->common->fl)); 4249566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 4259566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 426218db3c1SStefano Zampini #if defined(PETSC_USE_SUITESPARSE_GPU) 4279566063dSJacob 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)); 428218db3c1SStefano Zampini #endif 429641875f9SMatthew G Knepley 430641875f9SMatthew G Knepley F->ops->solve = MatSolve_CHOLMOD; 431641875f9SMatthew G Knepley F->ops->solvetranspose = MatSolve_CHOLMOD; 432218db3c1SStefano Zampini F->ops->matsolve = MatMatSolve_CHOLMOD; 433218db3c1SStefano Zampini F->ops->matsolvetranspose = MatMatSolve_CHOLMOD; 4343ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 435641875f9SMatthew G Knepley } 436641875f9SMatthew G Knepley 437d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F, Mat A, IS perm, const MatFactorInfo *info) 438d71ae5a4SJacob Faibussowitsch { 4396b8f6f9dSBarry Smith Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 440d0609cedSBarry Smith int err; 441641875f9SMatthew G Knepley cholmod_sparse cholA; 442218db3c1SStefano Zampini PetscBool aijalloc, valloc; 443641875f9SMatthew G Knepley PetscInt *fset = 0; 444641875f9SMatthew G Knepley size_t fsize = 0; 445641875f9SMatthew G Knepley 446641875f9SMatthew G Knepley PetscFunctionBegin; 44726cc229bSBarry Smith /* Set options to F */ 44826cc229bSBarry Smith PetscCall(CholmodSetOptions(F)); 44926cc229bSBarry Smith 4509566063dSJacob Faibussowitsch PetscCall((*chol->Wrap)(A, PETSC_FALSE, &cholA, &aijalloc, &valloc)); 451641875f9SMatthew G Knepley static_F = F; 452641875f9SMatthew G Knepley if (chol->factor) { 453d0609cedSBarry Smith err = !cholmod_X_resymbol(&cholA, fset, fsize, (int)chol->pack, chol->factor, chol->common); 454d0609cedSBarry Smith PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed with status %d", chol->common->status); 455641875f9SMatthew G Knepley } else if (perm) { 456641875f9SMatthew G Knepley const PetscInt *ip; 4579566063dSJacob Faibussowitsch PetscCall(ISGetIndices(perm, &ip)); 458641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze_p(&cholA, (PetscInt *)ip, fset, fsize, chol->common); 45928b400f6SJacob Faibussowitsch PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using PETSc ordering with status %d", chol->common->status); 4609566063dSJacob Faibussowitsch PetscCall(ISRestoreIndices(perm, &ip)); 461641875f9SMatthew G Knepley } else { 462641875f9SMatthew G Knepley chol->factor = cholmod_X_analyze(&cholA, chol->common); 46328b400f6SJacob Faibussowitsch PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status); 464641875f9SMatthew G Knepley } 465641875f9SMatthew G Knepley 4669566063dSJacob Faibussowitsch if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i)); 4679566063dSJacob Faibussowitsch if (valloc) PetscCall(PetscFree(cholA.x)); 468641875f9SMatthew G Knepley 469641875f9SMatthew G Knepley F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD; 4703ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 471641875f9SMatthew G Knepley } 472641875f9SMatthew G Knepley 473d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A, MatSolverType *type) 474d71ae5a4SJacob Faibussowitsch { 475641875f9SMatthew G Knepley PetscFunctionBegin; 476641875f9SMatthew G Knepley *type = MATSOLVERCHOLMOD; 4773ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 478641875f9SMatthew G Knepley } 479641875f9SMatthew G Knepley 480d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F, MatInfoType flag, MatInfo *info) 481d71ae5a4SJacob Faibussowitsch { 482218db3c1SStefano Zampini Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data; 483218db3c1SStefano Zampini 484218db3c1SStefano Zampini PetscFunctionBegin; 485218db3c1SStefano Zampini info->block_size = 1.0; 486218db3c1SStefano Zampini info->nz_allocated = chol->common->lnz; 487218db3c1SStefano Zampini info->nz_used = chol->common->lnz; 488218db3c1SStefano Zampini info->nz_unneeded = 0.0; 489218db3c1SStefano Zampini info->assemblies = 0.0; 490218db3c1SStefano Zampini info->mallocs = 0.0; 491218db3c1SStefano Zampini info->memory = chol->common->memory_inuse; 492218db3c1SStefano Zampini info->fill_ratio_given = 0; 493218db3c1SStefano Zampini info->fill_ratio_needed = 0; 494218db3c1SStefano Zampini info->factor_mallocs = chol->common->malloc_count; 4953ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 496218db3c1SStefano Zampini } 497218db3c1SStefano Zampini 498641875f9SMatthew G Knepley /*MC 499ee300463SSatish Balay MATSOLVERCHOLMOD 500ee300463SSatish Balay 501ee300463SSatish Balay A matrix type providing direct solvers (Cholesky) for sequential matrices 502641875f9SMatthew G Knepley via the external package CHOLMOD. 503641875f9SMatthew G Knepley 5042ef1f0ffSBarry Smith Use `./configure --download-suitesparse` to install PETSc to use CHOLMOD 505c2b89b5dSBarry Smith 5062ef1f0ffSBarry Smith Use `-pc_type cholesky` `-pc_factor_mat_solver_type cholmod` to use this direct solver 507641875f9SMatthew G Knepley 5082ef1f0ffSBarry Smith Consult CHOLMOD documentation for more information about the common parameters 509641875f9SMatthew G Knepley which correspond to the options database keys below. 510641875f9SMatthew G Knepley 511641875f9SMatthew G Knepley Options Database Keys: 512e08999f5SMatthew G Knepley + -mat_cholmod_dbound <0> - Minimum absolute value of diagonal entries of D (None) 513e08999f5SMatthew G Knepley . -mat_cholmod_grow0 <1.2> - Global growth ratio when factors are modified (None) 514e08999f5SMatthew G Knepley . -mat_cholmod_grow1 <1.2> - Column growth ratio when factors are modified (None) 515e08999f5SMatthew G Knepley . -mat_cholmod_grow2 <5> - Affine column growth constant when factors are modified (None) 516e08999f5SMatthew G Knepley . -mat_cholmod_maxrank <8> - Max rank of update, larger values are faster but use more memory [2,4,8] (None) 5172ef1f0ffSBarry Smith . -mat_cholmod_factor <AUTO> - (choose one of) `SIMPLICIAL`, `AUTO`, `SUPERNODAL` 518e08999f5SMatthew G Knepley . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None) 519e08999f5SMatthew G Knepley . -mat_cholmod_final_asis <TRUE> - Leave factors "as is" (None) 520e08999f5SMatthew G Knepley . -mat_cholmod_final_pack <TRUE> - Pack the columns when finished (use FALSE if the factors will be updated later) (None) 521e08999f5SMatthew G Knepley . -mat_cholmod_zrelax <0.8> - 3 real supernodal relaxed amalgamation parameters (None) 522e08999f5SMatthew G Knepley . -mat_cholmod_nrelax <4> - 3 size_t supernodal relaxed amalgamation parameters (None) 523e08999f5SMatthew G Knepley . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None) 5242c7c0729SBarry Smith . -mat_cholmod_print <3> - Verbosity level (None) 525ee300463SSatish Balay - -mat_ordering_type internal - Use the ordering provided by Cholmod 526641875f9SMatthew G Knepley 527641875f9SMatthew G Knepley Level: beginner 528641875f9SMatthew G Knepley 5292ef1f0ffSBarry Smith Note: 5302ef1f0ffSBarry Smith CHOLMOD is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html 531a364b7d2SBarry Smith 5322ef1f0ffSBarry Smith .seealso: [](chapter_matrices), `Mat`, `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType` 533641875f9SMatthew G Knepley M*/ 534b2573a8aSBarry Smith 535d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A, MatFactorType ftype, Mat *F) 536d71ae5a4SJacob Faibussowitsch { 537641875f9SMatthew G Knepley Mat B; 538641875f9SMatthew G Knepley Mat_CHOLMOD *chol; 539641875f9SMatthew G Knepley PetscInt m = A->rmap->n, n = A->cmap->n, bs; 540641875f9SMatthew G Knepley 541641875f9SMatthew G Knepley PetscFunctionBegin; 5429566063dSJacob Faibussowitsch PetscCall(MatGetBlockSize(A, &bs)); 54308401ef6SPierre Jolivet PetscCheck(bs == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "CHOLMOD only supports block size=1, given %" PetscInt_FMT, bs); 544218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX) 545b94d7dedSBarry Smith PetscCheck(A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only for Hermitian matrices"); 546218db3c1SStefano Zampini #endif 547641875f9SMatthew G Knepley /* Create the factorization matrix F */ 5489566063dSJacob Faibussowitsch PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B)); 5499566063dSJacob Faibussowitsch PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n)); 5509566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name)); 5519566063dSJacob Faibussowitsch PetscCall(MatSetUp(B)); 5524dfa11a4SJacob Faibussowitsch PetscCall(PetscNew(&chol)); 55326fbe8dcSKarl Rupp 554641875f9SMatthew G Knepley chol->Wrap = MatWrapCholmod_seqsbaij; 5556b8f6f9dSBarry Smith B->data = chol; 556641875f9SMatthew G Knepley 557218db3c1SStefano Zampini B->ops->getinfo = MatGetInfo_CHOLMOD; 558641875f9SMatthew G Knepley B->ops->view = MatView_CHOLMOD; 559641875f9SMatthew G Knepley B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD; 560641875f9SMatthew G Knepley B->ops->destroy = MatDestroy_CHOLMOD; 5619566063dSJacob Faibussowitsch PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqsbaij_cholmod)); 562641875f9SMatthew G Knepley B->factortype = MAT_FACTOR_CHOLESKY; 563218db3c1SStefano Zampini B->assembled = PETSC_TRUE; 564641875f9SMatthew G Knepley B->preallocated = PETSC_TRUE; 565641875f9SMatthew G Knepley 5669566063dSJacob Faibussowitsch PetscCall(CholmodStart(B)); 56700c67f3bSHong Zhang 5689566063dSJacob Faibussowitsch PetscCall(PetscFree(B->solvertype)); 5699566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype)); 570f73b0415SBarry Smith B->canuseordering = PETSC_TRUE; 5719566063dSJacob Faibussowitsch PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY])); 572641875f9SMatthew G Knepley *F = B; 5733ba16761SJacob Faibussowitsch PetscFunctionReturn(PETSC_SUCCESS); 574641875f9SMatthew G Knepley } 575