xref: /petsc/src/mat/impls/sbaij/seq/cholmod/sbaijcholmod.c (revision 48a46eb9bd028bec07ec0f396b1a3abb43f14558)
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
89e475b0dSSatish Balay    that one cannot use 64BIT_INDICES on 32bit machines [as Suitesparse_long is 32bit 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 
219371c9d4SSatish Balay static void CholmodErrorHandler(int status, const char *file, int line, const char *message) {
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 
629371c9d4SSatish Balay static PetscErrorCode CholmodSetOptions(Mat F) {
6326cc229bSBarry Smith   Mat_CHOLMOD    *chol = (Mat_CHOLMOD *)F->data;
6426cc229bSBarry Smith   cholmod_common *c    = chol->common;
6526cc229bSBarry Smith   PetscBool       flg;
6626cc229bSBarry Smith 
6726cc229bSBarry Smith   PetscFunctionBegin;
68d0609cedSBarry Smith   PetscOptionsBegin(PetscObjectComm((PetscObject)F), ((PetscObject)F)->prefix, "CHOLMOD Options", "Mat");
6954b3d318SStefano Zampini   CHOLMOD_OPTION_INT(nmethods, "Number of different ordering methods to try");
7026fbe8dcSKarl Rupp 
71b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU)
72b9eaa5e8SBarry Smith   c->useGPU = 1;
73b9eaa5e8SBarry Smith   CHOLMOD_OPTION_INT(useGPU, "Use GPU for BLAS 1, otherwise 0");
7454b3d318SStefano Zampini   CHOLMOD_OPTION_SIZE_T(maxGpuMemBytes, "Maximum memory to allocate on the GPU");
7554b3d318SStefano Zampini   CHOLMOD_OPTION_DOUBLE(maxGpuMemFraction, "Fraction of available GPU memory to allocate");
76b9eaa5e8SBarry Smith #endif
77b9eaa5e8SBarry Smith 
7854b3d318SStefano Zampini   /* CHOLMOD handles first-time packing and refactor-packing separately, but we usually want them to be the same. */
7954b3d318SStefano Zampini   chol->pack = (PetscBool)c->final_pack;
809566063dSJacob Faibussowitsch   PetscCall(PetscOptionsBool("-mat_cholmod_pack", "Pack factors after factorization [disable for frequent repeat factorization]", "None", chol->pack, &chol->pack, NULL));
81641875f9SMatthew G Knepley   c->final_pack = (int)chol->pack;
82641875f9SMatthew G Knepley 
83641875f9SMatthew G Knepley   CHOLMOD_OPTION_DOUBLE(dbound, "Minimum absolute value of diagonal entries of D");
84641875f9SMatthew G Knepley   CHOLMOD_OPTION_DOUBLE(grow0, "Global growth ratio when factors are modified");
85641875f9SMatthew G Knepley   CHOLMOD_OPTION_DOUBLE(grow1, "Column growth ratio when factors are modified");
86641875f9SMatthew G Knepley   CHOLMOD_OPTION_SIZE_T(grow2, "Affine column growth constant when factors are modified");
87641875f9SMatthew G Knepley   CHOLMOD_OPTION_SIZE_T(maxrank, "Max rank of update, larger values are faster but use more memory [2,4,8]");
88641875f9SMatthew G Knepley   {
89641875f9SMatthew G Knepley     static const char *const list[] = {"SIMPLICIAL", "AUTO", "SUPERNODAL", "MatCholmodFactorType", "MAT_CHOLMOD_FACTOR_", 0};
909566063dSJacob Faibussowitsch     PetscCall(PetscOptionsEnum("-mat_cholmod_factor", "Factorization method", "None", list, (PetscEnum)c->supernodal, (PetscEnum *)&c->supernodal, NULL));
91641875f9SMatthew G Knepley   }
92641875f9SMatthew G Knepley   if (c->supernodal) CHOLMOD_OPTION_DOUBLE(supernodal_switch, "flop/nnz_L threshold for switching to supernodal factorization");
93b9eaa5e8SBarry Smith   CHOLMOD_OPTION_BOOL(final_asis, "Leave factors \"as is\"");
94b9eaa5e8SBarry Smith   CHOLMOD_OPTION_BOOL(final_pack, "Pack the columns when finished (use FALSE if the factors will be updated later)");
95641875f9SMatthew G Knepley   if (!c->final_asis) {
96b9eaa5e8SBarry Smith     CHOLMOD_OPTION_BOOL(final_super, "Leave supernodal factors instead of converting to simplicial");
97b9eaa5e8SBarry Smith     CHOLMOD_OPTION_BOOL(final_ll, "Turn LDL' factorization into LL'");
98b9eaa5e8SBarry Smith     CHOLMOD_OPTION_BOOL(final_monotonic, "Ensure columns are monotonic when done");
99b9eaa5e8SBarry Smith     CHOLMOD_OPTION_BOOL(final_resymbol, "Remove numerically zero values resulting from relaxed supernodal amalgamation");
100641875f9SMatthew G Knepley   }
101641875f9SMatthew G Knepley   {
102641875f9SMatthew G Knepley     PetscReal tmp[] = {(PetscReal)c->zrelax[0], (PetscReal)c->zrelax[1], (PetscReal)c->zrelax[2]};
103641875f9SMatthew G Knepley     PetscInt  n     = 3;
1049566063dSJacob Faibussowitsch     PetscCall(PetscOptionsRealArray("-mat_cholmod_zrelax", "3 real supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg));
105aed4548fSBarry Smith     PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_zrelax");
1069371c9d4SSatish Balay     if (flg)
1079371c9d4SSatish Balay       while (n--) c->zrelax[n] = (double)tmp[n];
108641875f9SMatthew G Knepley   }
109641875f9SMatthew G Knepley   {
110641875f9SMatthew G Knepley     PetscInt n, tmp[] = {(PetscInt)c->nrelax[0], (PetscInt)c->nrelax[1], (PetscInt)c->nrelax[2]};
1119566063dSJacob Faibussowitsch     PetscCall(PetscOptionsIntArray("-mat_cholmod_nrelax", "3 size_t supernodal relaxed amalgamation parameters", "None", tmp, &n, &flg));
112aed4548fSBarry Smith     PetscCheck(!flg || n == 3, PetscObjectComm((PetscObject)F), PETSC_ERR_ARG_OUTOFRANGE, "must provide exactly 3 parameters to -mat_cholmod_nrelax");
1139371c9d4SSatish Balay     if (flg)
1149371c9d4SSatish Balay       while (n--) c->nrelax[n] = (size_t)tmp[n];
115641875f9SMatthew G Knepley   }
116b9eaa5e8SBarry Smith   CHOLMOD_OPTION_BOOL(prefer_upper, "Work with upper triangular form [faster when using fill-reducing ordering, slower in natural ordering]");
117b9eaa5e8SBarry Smith   CHOLMOD_OPTION_BOOL(default_nesdis, "Use NESDIS instead of METIS for nested dissection");
118641875f9SMatthew G Knepley   CHOLMOD_OPTION_INT(print, "Verbosity level");
119d0609cedSBarry Smith   PetscOptionsEnd();
120641875f9SMatthew G Knepley   PetscFunctionReturn(0);
121641875f9SMatthew G Knepley }
122641875f9SMatthew G Knepley 
1239371c9d4SSatish Balay PetscErrorCode CholmodStart(Mat F) {
12426cc229bSBarry Smith   Mat_CHOLMOD    *chol = (Mat_CHOLMOD *)F->data;
12526cc229bSBarry Smith   cholmod_common *c;
12626cc229bSBarry Smith 
12726cc229bSBarry Smith   PetscFunctionBegin;
12826cc229bSBarry Smith   if (chol->common) PetscFunctionReturn(0);
12926cc229bSBarry Smith   PetscCall(PetscMalloc1(1, &chol->common));
13026cc229bSBarry Smith   PetscCall(!cholmod_X_start(chol->common));
13126cc229bSBarry Smith 
13226cc229bSBarry Smith   c                = chol->common;
13326cc229bSBarry Smith   c->error_handler = CholmodErrorHandler;
13426cc229bSBarry Smith   PetscFunctionReturn(0);
13526cc229bSBarry Smith }
13626cc229bSBarry Smith 
1379371c9d4SSatish Balay static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc) {
138641875f9SMatthew G Knepley   Mat_SeqSBAIJ *sbaij    = (Mat_SeqSBAIJ *)A->data;
139218db3c1SStefano Zampini   PetscBool     vallocin = PETSC_FALSE;
140641875f9SMatthew G Knepley 
141641875f9SMatthew G Knepley   PetscFunctionBegin;
1429566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(C, sizeof(*C)));
143641875f9SMatthew 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 */
144641875f9SMatthew G Knepley   C->nrow  = (size_t)A->cmap->n;
145641875f9SMatthew G Knepley   C->ncol  = (size_t)A->rmap->n;
146641875f9SMatthew G Knepley   C->nzmax = (size_t)sbaij->maxnz;
147641875f9SMatthew G Knepley   C->p     = sbaij->i;
148641875f9SMatthew G Knepley   C->i     = sbaij->j;
149218db3c1SStefano Zampini   if (values) {
150218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
151218db3c1SStefano Zampini     /* we need to pass CHOLMOD the conjugate matrix */
152218db3c1SStefano Zampini     PetscScalar *v;
153218db3c1SStefano Zampini     PetscInt     i;
154218db3c1SStefano Zampini 
1559566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sbaij->maxnz, &v));
156218db3c1SStefano Zampini     for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]);
157218db3c1SStefano Zampini     C->x     = v;
158218db3c1SStefano Zampini     vallocin = PETSC_TRUE;
159218db3c1SStefano Zampini #else
160641875f9SMatthew G Knepley     C->x = sbaij->a;
161218db3c1SStefano Zampini #endif
162218db3c1SStefano Zampini   }
163641875f9SMatthew G Knepley   C->stype  = -1;
164641875f9SMatthew G Knepley   C->itype  = CHOLMOD_INT_TYPE;
165218db3c1SStefano Zampini   C->xtype  = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
166641875f9SMatthew G Knepley   C->dtype  = CHOLMOD_DOUBLE;
167641875f9SMatthew G Knepley   C->sorted = 1;
168641875f9SMatthew G Knepley   C->packed = 1;
169641875f9SMatthew G Knepley   *aijalloc = PETSC_FALSE;
170218db3c1SStefano Zampini   *valloc   = vallocin;
171641875f9SMatthew G Knepley   PetscFunctionReturn(0);
172641875f9SMatthew G Knepley }
173641875f9SMatthew G Knepley 
174218db3c1SStefano Zampini #define GET_ARRAY_READ  0
175218db3c1SStefano Zampini #define GET_ARRAY_WRITE 1
176218db3c1SStefano Zampini 
1779371c9d4SSatish Balay PetscErrorCode VecWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) {
178218db3c1SStefano Zampini   PetscScalar *x;
179641875f9SMatthew G Knepley   PetscInt     n;
180641875f9SMatthew G Knepley 
181641875f9SMatthew G Knepley   PetscFunctionBegin;
1829566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(Y, sizeof(*Y)));
183218db3c1SStefano Zampini   switch (rw) {
1849371c9d4SSatish Balay   case GET_ARRAY_READ: PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x)); break;
1859371c9d4SSatish Balay   case GET_ARRAY_WRITE: PetscCall(VecGetArrayWrite(X, &x)); break;
1869371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break;
187218db3c1SStefano Zampini   }
1889566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &n));
18926fbe8dcSKarl Rupp 
190218db3c1SStefano Zampini   Y->x     = x;
191641875f9SMatthew G Knepley   Y->nrow  = n;
192641875f9SMatthew G Knepley   Y->ncol  = 1;
193641875f9SMatthew G Knepley   Y->nzmax = n;
194641875f9SMatthew G Knepley   Y->d     = n;
195641875f9SMatthew G Knepley   Y->xtype = CHOLMOD_SCALAR_TYPE;
196641875f9SMatthew G Knepley   Y->dtype = CHOLMOD_DOUBLE;
197641875f9SMatthew G Knepley   PetscFunctionReturn(0);
198641875f9SMatthew G Knepley }
199641875f9SMatthew G Knepley 
2009371c9d4SSatish Balay PetscErrorCode VecUnWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y) {
201d9ca1df4SBarry Smith   PetscFunctionBegin;
202218db3c1SStefano Zampini   switch (rw) {
2039371c9d4SSatish Balay   case GET_ARRAY_READ: PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&Y->x)); break;
2049371c9d4SSatish Balay   case GET_ARRAY_WRITE: PetscCall(VecRestoreArrayWrite(X, (PetscScalar **)&Y->x)); break;
2059371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break;
206218db3c1SStefano Zampini   }
207218db3c1SStefano Zampini   PetscFunctionReturn(0);
208218db3c1SStefano Zampini }
209218db3c1SStefano Zampini 
2109371c9d4SSatish Balay PetscErrorCode MatDenseWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) {
211218db3c1SStefano Zampini   PetscScalar *x;
212218db3c1SStefano Zampini   PetscInt     m, n, lda;
213218db3c1SStefano Zampini 
214218db3c1SStefano Zampini   PetscFunctionBegin;
2159566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(Y, sizeof(*Y)));
216218db3c1SStefano Zampini   switch (rw) {
2179371c9d4SSatish Balay   case GET_ARRAY_READ: PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&x)); break;
2189371c9d4SSatish Balay   case GET_ARRAY_WRITE: PetscCall(MatDenseGetArrayWrite(X, &x)); break;
2199371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break;
220218db3c1SStefano Zampini   }
2219566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
2229566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, &n));
223218db3c1SStefano Zampini 
224218db3c1SStefano Zampini   Y->x     = x;
225218db3c1SStefano Zampini   Y->nrow  = m;
226218db3c1SStefano Zampini   Y->ncol  = n;
227218db3c1SStefano Zampini   Y->nzmax = lda * n;
228218db3c1SStefano Zampini   Y->d     = lda;
229218db3c1SStefano Zampini   Y->xtype = CHOLMOD_SCALAR_TYPE;
230218db3c1SStefano Zampini   Y->dtype = CHOLMOD_DOUBLE;
231218db3c1SStefano Zampini   PetscFunctionReturn(0);
232218db3c1SStefano Zampini }
233218db3c1SStefano Zampini 
2349371c9d4SSatish Balay PetscErrorCode MatDenseUnWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y) {
235218db3c1SStefano Zampini   PetscFunctionBegin;
236218db3c1SStefano Zampini   switch (rw) {
2379371c9d4SSatish Balay   case GET_ARRAY_READ: PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&Y->x)); break;
238218db3c1SStefano Zampini   case GET_ARRAY_WRITE:
239218db3c1SStefano Zampini     /* we don't have MatDenseRestoreArrayWrite */
2409566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X, (PetscScalar **)&Y->x));
241218db3c1SStefano Zampini     break;
2429371c9d4SSatish Balay   default: SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw); break;
243218db3c1SStefano Zampini   }
244d9ca1df4SBarry Smith   PetscFunctionReturn(0);
245d9ca1df4SBarry Smith }
246d9ca1df4SBarry Smith 
2479371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F) {
2486b8f6f9dSBarry Smith   Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
249641875f9SMatthew G Knepley 
250641875f9SMatthew G Knepley   PetscFunctionBegin;
251*48a46eb9SPierre Jolivet   if (chol->spqrfact) PetscCall(!SuiteSparseQR_C_free(&chol->spqrfact, chol->common));
252*48a46eb9SPierre Jolivet   if (chol->factor) PetscCall(!cholmod_X_free_factor(&chol->factor, chol->common));
253a2fc1e05SToby Isaac   if (chol->common->itype == CHOLMOD_INT) {
2549566063dSJacob Faibussowitsch     PetscCall(!cholmod_finish(chol->common));
255a2fc1e05SToby Isaac   } else {
2569566063dSJacob Faibussowitsch     PetscCall(!cholmod_l_finish(chol->common));
257a2fc1e05SToby Isaac   }
2589566063dSJacob Faibussowitsch   PetscCall(PetscFree(chol->common));
2599566063dSJacob Faibussowitsch   PetscCall(PetscFree(chol->matrix));
2609566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", NULL));
2612e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorSymbolic_C", NULL));
2622e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", NULL));
2639566063dSJacob Faibussowitsch   PetscCall(PetscFree(F->data));
264641875f9SMatthew G Knepley   PetscFunctionReturn(0);
265641875f9SMatthew G Knepley }
266641875f9SMatthew G Knepley 
267641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat, Vec, Vec);
268218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat, Mat, Mat);
269641875f9SMatthew G Knepley 
270fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/
271641875f9SMatthew G Knepley 
2729371c9d4SSatish Balay static PetscErrorCode MatView_Info_CHOLMOD(Mat F, PetscViewer viewer) {
2736b8f6f9dSBarry Smith   Mat_CHOLMOD          *chol = (Mat_CHOLMOD *)F->data;
274641875f9SMatthew G Knepley   const cholmod_common *c    = chol->common;
275641875f9SMatthew G Knepley   PetscInt              i;
276641875f9SMatthew G Knepley 
277641875f9SMatthew G Knepley   PetscFunctionBegin;
278641875f9SMatthew G Knepley   if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(0);
2799566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "CHOLMOD run parameters:\n"));
2809566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
2819566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Pack factors after symbolic factorization: %s\n", chol->pack ? "TRUE" : "FALSE"));
2829566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.dbound            %g  (Smallest absolute value of diagonal entries of D)\n", c->dbound));
2839566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow0             %g\n", c->grow0));
2849566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow1             %g\n", c->grow1));
2859566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow2             %u\n", (unsigned)c->grow2));
2869566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.maxrank           %u\n", (unsigned)c->maxrank));
2879566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal_switch %g\n", c->supernodal_switch));
2889566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal        %d\n", c->supernodal));
2899566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_asis        %d\n", c->final_asis));
2909566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_super       %d\n", c->final_super));
2919566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_ll          %d\n", c->final_ll));
2929566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_pack        %d\n", c->final_pack));
2939566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_monotonic   %d\n", c->final_monotonic));
2949566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_resymbol    %d\n", c->final_resymbol));
2959566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.zrelax            [%g,%g,%g]\n", c->zrelax[0], c->zrelax[1], c->zrelax[2]));
2969566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrelax            [%u,%u,%u]\n", (unsigned)c->nrelax[0], (unsigned)c->nrelax[1], (unsigned)c->nrelax[2]));
2979566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.prefer_upper      %d\n", c->prefer_upper));
2989566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.print             %d\n", c->print));
299641875f9SMatthew G Knepley   for (i = 0; i < c->nmethods; i++) {
3009566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering method %" PetscInt_FMT "%s:\n", i, i == c->selected ? " [SELECTED]" : ""));
3019371c9d4SSatish 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));
302641875f9SMatthew G Knepley   }
3039566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.postorder         %d\n", c->postorder));
3049566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.default_nesdis    %d (use NESDIS instead of METIS for nested dissection)\n", c->default_nesdis));
305641875f9SMatthew G Knepley   /* Statistics */
3069566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.fl                %g (flop count from most recent analysis)\n", c->fl));
3079566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.lnz               %g (fundamental nz in L)\n", c->lnz));
3089566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.anz               %g\n", c->anz));
3099566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.modfl             %g (flop count from most recent update)\n", c->modfl));
3109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.malloc_count      %g (number of live objects)\n", (double)c->malloc_count));
3119566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_usage      %g (peak memory usage in bytes)\n", (double)c->memory_usage));
3129566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_inuse      %g (current memory usage in bytes)\n", (double)c->memory_inuse));
3139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_col      %g (number of column reallocations)\n", c->nrealloc_col));
3149566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_factor   %g (number of factor reallocations due to column reallocations)\n", c->nrealloc_factor));
3159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.ndbounds_hit      %g (number of times diagonal was modified by dbound)\n", c->ndbounds_hit));
3169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.rowfacfl          %g (number of flops in last call to cholmod_rowfac)\n", c->rowfacfl));
3179566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.aatfl             %g (number of flops to compute A(:,f)*A(:,f)')\n", c->aatfl));
318b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU)
3199566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.useGPU            %d\n", c->useGPU));
320b9eaa5e8SBarry Smith #endif
3219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
322641875f9SMatthew G Knepley   PetscFunctionReturn(0);
323641875f9SMatthew G Knepley }
324641875f9SMatthew G Knepley 
3259371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F, PetscViewer viewer) {
326ace3abfcSBarry Smith   PetscBool         iascii;
327641875f9SMatthew G Knepley   PetscViewerFormat format;
328641875f9SMatthew G Knepley 
329641875f9SMatthew G Knepley   PetscFunctionBegin;
3309566063dSJacob Faibussowitsch   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &iascii));
331641875f9SMatthew G Knepley   if (iascii) {
3329566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
333*48a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_CHOLMOD(F, viewer));
334641875f9SMatthew G Knepley   }
335641875f9SMatthew G Knepley   PetscFunctionReturn(0);
336641875f9SMatthew G Knepley }
337641875f9SMatthew G Knepley 
3389371c9d4SSatish Balay static PetscErrorCode MatSolve_CHOLMOD(Mat F, Vec B, Vec X) {
3396b8f6f9dSBarry Smith   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
340218db3c1SStefano Zampini   cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL;
341641875f9SMatthew G Knepley 
342641875f9SMatthew G Knepley   PetscFunctionBegin;
343641875f9SMatthew G Knepley   static_F = F;
3449566063dSJacob Faibussowitsch   PetscCall(VecWrapCholmod(B, GET_ARRAY_READ, &cholB));
3459566063dSJacob Faibussowitsch   PetscCall(VecWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
346218db3c1SStefano Zampini   X_handle = &cholX;
3479566063dSJacob Faibussowitsch   PetscCall(!cholmod_X_solve2(CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common));
3489566063dSJacob Faibussowitsch   PetscCall(!cholmod_X_free_dense(&Y_handle, chol->common));
3499566063dSJacob Faibussowitsch   PetscCall(!cholmod_X_free_dense(&E_handle, chol->common));
3509566063dSJacob Faibussowitsch   PetscCall(VecUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
3519566063dSJacob Faibussowitsch   PetscCall(VecUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
3529566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(4.0 * chol->common->lnz));
353218db3c1SStefano Zampini   PetscFunctionReturn(0);
354218db3c1SStefano Zampini }
355218db3c1SStefano Zampini 
3569371c9d4SSatish Balay static PetscErrorCode MatMatSolve_CHOLMOD(Mat F, Mat B, Mat X) {
357218db3c1SStefano Zampini   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
358218db3c1SStefano Zampini   cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL;
359218db3c1SStefano Zampini 
360218db3c1SStefano Zampini   PetscFunctionBegin;
361218db3c1SStefano Zampini   static_F = F;
3629566063dSJacob Faibussowitsch   PetscCall(MatDenseWrapCholmod(B, GET_ARRAY_READ, &cholB));
3639566063dSJacob Faibussowitsch   PetscCall(MatDenseWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
364218db3c1SStefano Zampini   X_handle = &cholX;
3659566063dSJacob Faibussowitsch   PetscCall(!cholmod_X_solve2(CHOLMOD_A, chol->factor, &cholB, NULL, &X_handle, NULL, &Y_handle, &E_handle, chol->common));
3669566063dSJacob Faibussowitsch   PetscCall(!cholmod_X_free_dense(&Y_handle, chol->common));
3679566063dSJacob Faibussowitsch   PetscCall(!cholmod_X_free_dense(&E_handle, chol->common));
3689566063dSJacob Faibussowitsch   PetscCall(MatDenseUnWrapCholmod(B, GET_ARRAY_READ, &cholB));
3699566063dSJacob Faibussowitsch   PetscCall(MatDenseUnWrapCholmod(X, GET_ARRAY_WRITE, &cholX));
3709566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(4.0 * B->cmap->n * chol->common->lnz));
371641875f9SMatthew G Knepley   PetscFunctionReturn(0);
372641875f9SMatthew G Knepley }
373641875f9SMatthew G Knepley 
3749371c9d4SSatish Balay static PetscErrorCode MatCholeskyFactorNumeric_CHOLMOD(Mat F, Mat A, const MatFactorInfo *info) {
3756b8f6f9dSBarry Smith   Mat_CHOLMOD   *chol = (Mat_CHOLMOD *)F->data;
376641875f9SMatthew G Knepley   cholmod_sparse cholA;
377218db3c1SStefano Zampini   PetscBool      aijalloc, valloc;
378d0609cedSBarry Smith   int            err;
379641875f9SMatthew G Knepley 
380641875f9SMatthew G Knepley   PetscFunctionBegin;
3819566063dSJacob Faibussowitsch   PetscCall((*chol->Wrap)(A, PETSC_TRUE, &cholA, &aijalloc, &valloc));
382641875f9SMatthew G Knepley   static_F = F;
383d0609cedSBarry Smith   err      = !cholmod_X_factorize(&cholA, chol->factor, chol->common);
384d0609cedSBarry Smith   PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD factorization failed with status %d", chol->common->status);
38508401ef6SPierre 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);
386641875f9SMatthew G Knepley 
3879566063dSJacob Faibussowitsch   PetscCall(PetscLogFlops(chol->common->fl));
3889566063dSJacob Faibussowitsch   if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
3899566063dSJacob Faibussowitsch   if (valloc) PetscCall(PetscFree(cholA.x));
390218db3c1SStefano Zampini #if defined(PETSC_USE_SUITESPARSE_GPU)
3919566063dSJacob 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));
392218db3c1SStefano Zampini #endif
393641875f9SMatthew G Knepley 
394641875f9SMatthew G Knepley   F->ops->solve             = MatSolve_CHOLMOD;
395641875f9SMatthew G Knepley   F->ops->solvetranspose    = MatSolve_CHOLMOD;
396218db3c1SStefano Zampini   F->ops->matsolve          = MatMatSolve_CHOLMOD;
397218db3c1SStefano Zampini   F->ops->matsolvetranspose = MatMatSolve_CHOLMOD;
398641875f9SMatthew G Knepley   PetscFunctionReturn(0);
399641875f9SMatthew G Knepley }
400641875f9SMatthew G Knepley 
4019371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatCholeskyFactorSymbolic_CHOLMOD(Mat F, Mat A, IS perm, const MatFactorInfo *info) {
4026b8f6f9dSBarry Smith   Mat_CHOLMOD   *chol = (Mat_CHOLMOD *)F->data;
403d0609cedSBarry Smith   int            err;
404641875f9SMatthew G Knepley   cholmod_sparse cholA;
405218db3c1SStefano Zampini   PetscBool      aijalloc, valloc;
406641875f9SMatthew G Knepley   PetscInt      *fset  = 0;
407641875f9SMatthew G Knepley   size_t         fsize = 0;
408641875f9SMatthew G Knepley 
409641875f9SMatthew G Knepley   PetscFunctionBegin;
41026cc229bSBarry Smith   /* Set options to F */
41126cc229bSBarry Smith   PetscCall(CholmodSetOptions(F));
41226cc229bSBarry Smith 
4139566063dSJacob Faibussowitsch   PetscCall((*chol->Wrap)(A, PETSC_FALSE, &cholA, &aijalloc, &valloc));
414641875f9SMatthew G Knepley   static_F = F;
415641875f9SMatthew G Knepley   if (chol->factor) {
416d0609cedSBarry Smith     err = !cholmod_X_resymbol(&cholA, fset, fsize, (int)chol->pack, chol->factor, chol->common);
417d0609cedSBarry Smith     PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed with status %d", chol->common->status);
418641875f9SMatthew G Knepley   } else if (perm) {
419641875f9SMatthew G Knepley     const PetscInt *ip;
4209566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(perm, &ip));
421641875f9SMatthew G Knepley     chol->factor = cholmod_X_analyze_p(&cholA, (PetscInt *)ip, fset, fsize, chol->common);
42228b400f6SJacob Faibussowitsch     PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using PETSc ordering with status %d", chol->common->status);
4239566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(perm, &ip));
424641875f9SMatthew G Knepley   } else {
425641875f9SMatthew G Knepley     chol->factor = cholmod_X_analyze(&cholA, chol->common);
42628b400f6SJacob Faibussowitsch     PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status);
427641875f9SMatthew G Knepley   }
428641875f9SMatthew G Knepley 
4299566063dSJacob Faibussowitsch   if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
4309566063dSJacob Faibussowitsch   if (valloc) PetscCall(PetscFree(cholA.x));
431641875f9SMatthew G Knepley 
432641875f9SMatthew G Knepley   F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD;
433641875f9SMatthew G Knepley   PetscFunctionReturn(0);
434641875f9SMatthew G Knepley }
435641875f9SMatthew G Knepley 
4369371c9d4SSatish Balay static PetscErrorCode MatFactorGetSolverType_seqsbaij_cholmod(Mat A, MatSolverType *type) {
437641875f9SMatthew G Knepley   PetscFunctionBegin;
438641875f9SMatthew G Knepley   *type = MATSOLVERCHOLMOD;
439641875f9SMatthew G Knepley   PetscFunctionReturn(0);
440641875f9SMatthew G Knepley }
441641875f9SMatthew G Knepley 
4429371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatGetInfo_CHOLMOD(Mat F, MatInfoType flag, MatInfo *info) {
443218db3c1SStefano Zampini   Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
444218db3c1SStefano Zampini 
445218db3c1SStefano Zampini   PetscFunctionBegin;
446218db3c1SStefano Zampini   info->block_size        = 1.0;
447218db3c1SStefano Zampini   info->nz_allocated      = chol->common->lnz;
448218db3c1SStefano Zampini   info->nz_used           = chol->common->lnz;
449218db3c1SStefano Zampini   info->nz_unneeded       = 0.0;
450218db3c1SStefano Zampini   info->assemblies        = 0.0;
451218db3c1SStefano Zampini   info->mallocs           = 0.0;
452218db3c1SStefano Zampini   info->memory            = chol->common->memory_inuse;
453218db3c1SStefano Zampini   info->fill_ratio_given  = 0;
454218db3c1SStefano Zampini   info->fill_ratio_needed = 0;
455218db3c1SStefano Zampini   info->factor_mallocs    = chol->common->malloc_count;
456218db3c1SStefano Zampini   PetscFunctionReturn(0);
457218db3c1SStefano Zampini }
458218db3c1SStefano Zampini 
459641875f9SMatthew G Knepley /*MC
460ee300463SSatish Balay   MATSOLVERCHOLMOD
461ee300463SSatish Balay 
462ee300463SSatish Balay   A matrix type providing direct solvers (Cholesky) for sequential matrices
463641875f9SMatthew G Knepley   via the external package CHOLMOD.
464641875f9SMatthew G Knepley 
465c2b89b5dSBarry Smith   Use ./configure --download-suitesparse to install PETSc to use CHOLMOD
466c2b89b5dSBarry Smith 
467f29f8b16SBarry Smith   Use -pc_type cholesky -pc_factor_mat_solver_type cholmod to use this direct solver
468641875f9SMatthew G Knepley 
469641875f9SMatthew G Knepley   Consult CHOLMOD documentation for more information about the Common parameters
470641875f9SMatthew G Knepley   which correspond to the options database keys below.
471641875f9SMatthew G Knepley 
472641875f9SMatthew G Knepley   Options Database Keys:
473e08999f5SMatthew G Knepley + -mat_cholmod_dbound <0>          - Minimum absolute value of diagonal entries of D (None)
474e08999f5SMatthew G Knepley . -mat_cholmod_grow0 <1.2>         - Global growth ratio when factors are modified (None)
475e08999f5SMatthew G Knepley . -mat_cholmod_grow1 <1.2>         - Column growth ratio when factors are modified (None)
476e08999f5SMatthew G Knepley . -mat_cholmod_grow2 <5>           - Affine column growth constant when factors are modified (None)
477e08999f5SMatthew G Knepley . -mat_cholmod_maxrank <8>         - Max rank of update, larger values are faster but use more memory [2,4,8] (None)
478e08999f5SMatthew G Knepley . -mat_cholmod_factor <AUTO>       - (choose one of) SIMPLICIAL AUTO SUPERNODAL
479e08999f5SMatthew G Knepley . -mat_cholmod_supernodal_switch <40> - flop/nnz_L threshold for switching to supernodal factorization (None)
480e08999f5SMatthew G Knepley . -mat_cholmod_final_asis <TRUE>   - Leave factors "as is" (None)
481e08999f5SMatthew G Knepley . -mat_cholmod_final_pack <TRUE>   - Pack the columns when finished (use FALSE if the factors will be updated later) (None)
482e08999f5SMatthew G Knepley . -mat_cholmod_zrelax <0.8>        - 3 real supernodal relaxed amalgamation parameters (None)
483e08999f5SMatthew G Knepley . -mat_cholmod_nrelax <4>          - 3 size_t supernodal relaxed amalgamation parameters (None)
484e08999f5SMatthew G Knepley . -mat_cholmod_prefer_upper <TRUE> - Work with upper triangular form (faster when using fill-reducing ordering, slower in natural ordering) (None)
4852c7c0729SBarry Smith . -mat_cholmod_print <3>           - Verbosity level (None)
486ee300463SSatish Balay - -mat_ordering_type internal      - Use the ordering provided by Cholmod
487641875f9SMatthew G Knepley 
488641875f9SMatthew G Knepley    Level: beginner
489641875f9SMatthew G Knepley 
490a364b7d2SBarry Smith    Note: CHOLMOD is part of SuiteSparse http://faculty.cse.tamu.edu/davis/suitesparse.html
491a364b7d2SBarry Smith 
492db781477SPatrick Sanan .seealso: `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType`
493641875f9SMatthew G Knepley M*/
494b2573a8aSBarry Smith 
4959371c9d4SSatish Balay PETSC_INTERN PetscErrorCode MatGetFactor_seqsbaij_cholmod(Mat A, MatFactorType ftype, Mat *F) {
496641875f9SMatthew G Knepley   Mat          B;
497641875f9SMatthew G Knepley   Mat_CHOLMOD *chol;
498641875f9SMatthew G Knepley   PetscInt     m = A->rmap->n, n = A->cmap->n, bs;
499641875f9SMatthew G Knepley 
500641875f9SMatthew G Knepley   PetscFunctionBegin;
5019566063dSJacob Faibussowitsch   PetscCall(MatGetBlockSize(A, &bs));
50208401ef6SPierre Jolivet   PetscCheck(bs == 1, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "CHOLMOD only supports block size=1, given %" PetscInt_FMT, bs);
503218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
504b94d7dedSBarry Smith   PetscCheck(A->hermitian == PETSC_BOOL3_TRUE, PetscObjectComm((PetscObject)A), PETSC_ERR_SUP, "Only for Hermitian matrices");
505218db3c1SStefano Zampini #endif
506641875f9SMatthew G Knepley   /* Create the factorization matrix F */
5079566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5089566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
5099566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name));
5109566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
5119566063dSJacob Faibussowitsch   PetscCall(PetscNewLog(B, &chol));
51226fbe8dcSKarl Rupp 
513641875f9SMatthew G Knepley   chol->Wrap = MatWrapCholmod_seqsbaij;
5146b8f6f9dSBarry Smith   B->data    = chol;
515641875f9SMatthew G Knepley 
516218db3c1SStefano Zampini   B->ops->getinfo                = MatGetInfo_CHOLMOD;
517641875f9SMatthew G Knepley   B->ops->view                   = MatView_CHOLMOD;
518641875f9SMatthew G Knepley   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
519641875f9SMatthew G Knepley   B->ops->destroy                = MatDestroy_CHOLMOD;
5209566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqsbaij_cholmod));
521641875f9SMatthew G Knepley   B->factortype   = MAT_FACTOR_CHOLESKY;
522218db3c1SStefano Zampini   B->assembled    = PETSC_TRUE;
523641875f9SMatthew G Knepley   B->preallocated = PETSC_TRUE;
524641875f9SMatthew G Knepley 
5259566063dSJacob Faibussowitsch   PetscCall(CholmodStart(B));
52600c67f3bSHong Zhang 
5279566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
5289566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype));
529f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
5309566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
531641875f9SMatthew G Knepley   *F = B;
532641875f9SMatthew G Knepley   PetscFunctionReturn(0);
533641875f9SMatthew G Knepley }
534