xref: /petsc/src/mat/impls/sbaij/seq/cholmod/sbaijcholmod.c (revision ae1ee55146a7ad071171b861759b85940c7e5c67)
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 
CholmodErrorHandler(int status,const char * file,int line,const char * message)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; \
5052ce0ab5SPierre Jolivet     PetscCall(PetscOptionsBoundedReal("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL, 0.0)); \
51641875f9SMatthew G Knepley     c->name = (size_t)tmp; \
52641875f9SMatthew G Knepley   } while (0)
5326fbe8dcSKarl Rupp 
549371c9d4SSatish Balay #define CHOLMOD_OPTION_BOOL(name, help) \
559371c9d4SSatish Balay   do { \
56ace3abfcSBarry Smith     PetscBool tmp = (PetscBool)!!c->name; \
579566063dSJacob Faibussowitsch     PetscCall(PetscOptionsBool("-mat_cholmod_" #name, help, "None", tmp, &tmp, NULL)); \
58641875f9SMatthew G Knepley     c->name = (int)tmp; \
59641875f9SMatthew G Knepley   } while (0)
60641875f9SMatthew G Knepley 
CholmodSetOptions(Mat F)61d71ae5a4SJacob Faibussowitsch static PetscErrorCode CholmodSetOptions(Mat F)
62d71ae5a4SJacob Faibussowitsch {
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();
1203ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
121641875f9SMatthew G Knepley }
122641875f9SMatthew G Knepley 
CholmodStart(Mat F)123d71ae5a4SJacob Faibussowitsch PetscErrorCode CholmodStart(Mat F)
124d71ae5a4SJacob Faibussowitsch {
12526cc229bSBarry Smith   Mat_CHOLMOD    *chol = (Mat_CHOLMOD *)F->data;
12626cc229bSBarry Smith   cholmod_common *c;
12726cc229bSBarry Smith 
12826cc229bSBarry Smith   PetscFunctionBegin;
1293ba16761SJacob Faibussowitsch   if (chol->common) PetscFunctionReturn(PETSC_SUCCESS);
13026cc229bSBarry Smith   PetscCall(PetscMalloc1(1, &chol->common));
1313ba16761SJacob Faibussowitsch   PetscCallExternal(!cholmod_X_start, chol->common);
13226cc229bSBarry Smith 
13326cc229bSBarry Smith   c                = chol->common;
13426cc229bSBarry Smith   c->error_handler = CholmodErrorHandler;
1353ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
13626cc229bSBarry Smith }
13726cc229bSBarry Smith 
MatWrapCholmod_seqsbaij(Mat A,PetscBool values,cholmod_sparse * C,PetscBool * aijalloc,PetscBool * valloc)138d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatWrapCholmod_seqsbaij(Mat A, PetscBool values, cholmod_sparse *C, PetscBool *aijalloc, PetscBool *valloc)
139d71ae5a4SJacob Faibussowitsch {
140641875f9SMatthew G Knepley   Mat_SeqSBAIJ *sbaij    = (Mat_SeqSBAIJ *)A->data;
141218db3c1SStefano Zampini   PetscBool     vallocin = PETSC_FALSE;
142641875f9SMatthew G Knepley 
143641875f9SMatthew G Knepley   PetscFunctionBegin;
1449566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(C, sizeof(*C)));
145641875f9SMatthew 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 */
146641875f9SMatthew G Knepley   C->nrow  = (size_t)A->cmap->n;
147641875f9SMatthew G Knepley   C->ncol  = (size_t)A->rmap->n;
148641875f9SMatthew G Knepley   C->nzmax = (size_t)sbaij->maxnz;
149641875f9SMatthew G Knepley   C->p     = sbaij->i;
150641875f9SMatthew G Knepley   C->i     = sbaij->j;
151218db3c1SStefano Zampini   if (values) {
152218db3c1SStefano Zampini #if defined(PETSC_USE_COMPLEX)
153218db3c1SStefano Zampini     /* we need to pass CHOLMOD the conjugate matrix */
154218db3c1SStefano Zampini     PetscScalar *v;
155218db3c1SStefano Zampini     PetscInt     i;
156218db3c1SStefano Zampini 
1579566063dSJacob Faibussowitsch     PetscCall(PetscMalloc1(sbaij->maxnz, &v));
158218db3c1SStefano Zampini     for (i = 0; i < sbaij->maxnz; i++) v[i] = PetscConj(sbaij->a[i]);
159218db3c1SStefano Zampini     C->x     = v;
160218db3c1SStefano Zampini     vallocin = PETSC_TRUE;
161218db3c1SStefano Zampini #else
162641875f9SMatthew G Knepley     C->x = sbaij->a;
163218db3c1SStefano Zampini #endif
164218db3c1SStefano Zampini   }
165641875f9SMatthew G Knepley   C->stype  = -1;
166641875f9SMatthew G Knepley   C->itype  = CHOLMOD_INT_TYPE;
167218db3c1SStefano Zampini   C->xtype  = values ? CHOLMOD_SCALAR_TYPE : CHOLMOD_PATTERN;
168641875f9SMatthew G Knepley   C->dtype  = CHOLMOD_DOUBLE;
169641875f9SMatthew G Knepley   C->sorted = 1;
170641875f9SMatthew G Knepley   C->packed = 1;
171641875f9SMatthew G Knepley   *aijalloc = PETSC_FALSE;
172218db3c1SStefano Zampini   *valloc   = vallocin;
1733ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
174641875f9SMatthew G Knepley }
175641875f9SMatthew G Knepley 
176218db3c1SStefano Zampini #define GET_ARRAY_READ  0
177218db3c1SStefano Zampini #define GET_ARRAY_WRITE 1
178218db3c1SStefano Zampini 
VecWrapCholmod(Vec X,PetscInt rw,cholmod_dense * Y)179d71ae5a4SJacob Faibussowitsch PetscErrorCode VecWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y)
180d71ae5a4SJacob Faibussowitsch {
181218db3c1SStefano Zampini   PetscScalar *x;
182641875f9SMatthew G Knepley   PetscInt     n;
183641875f9SMatthew G Knepley 
184641875f9SMatthew G Knepley   PetscFunctionBegin;
1859566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(Y, sizeof(*Y)));
186218db3c1SStefano Zampini   switch (rw) {
187d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_READ:
188d71ae5a4SJacob Faibussowitsch     PetscCall(VecGetArrayRead(X, (const PetscScalar **)&x));
189d71ae5a4SJacob Faibussowitsch     break;
190d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_WRITE:
191d71ae5a4SJacob Faibussowitsch     PetscCall(VecGetArrayWrite(X, &x));
192d71ae5a4SJacob Faibussowitsch     break;
193d71ae5a4SJacob Faibussowitsch   default:
194d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
195d71ae5a4SJacob Faibussowitsch     break;
196218db3c1SStefano Zampini   }
1979566063dSJacob Faibussowitsch   PetscCall(VecGetSize(X, &n));
19826fbe8dcSKarl Rupp 
199218db3c1SStefano Zampini   Y->x     = x;
200641875f9SMatthew G Knepley   Y->nrow  = n;
201641875f9SMatthew G Knepley   Y->ncol  = 1;
202641875f9SMatthew G Knepley   Y->nzmax = n;
203641875f9SMatthew G Knepley   Y->d     = n;
204641875f9SMatthew G Knepley   Y->xtype = CHOLMOD_SCALAR_TYPE;
205641875f9SMatthew G Knepley   Y->dtype = CHOLMOD_DOUBLE;
2063ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
207641875f9SMatthew G Knepley }
208641875f9SMatthew G Knepley 
VecUnWrapCholmod(Vec X,PetscInt rw,cholmod_dense * Y)209d71ae5a4SJacob Faibussowitsch PetscErrorCode VecUnWrapCholmod(Vec X, PetscInt rw, cholmod_dense *Y)
210d71ae5a4SJacob Faibussowitsch {
211d9ca1df4SBarry Smith   PetscFunctionBegin;
212218db3c1SStefano Zampini   switch (rw) {
213d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_READ:
214d71ae5a4SJacob Faibussowitsch     PetscCall(VecRestoreArrayRead(X, (const PetscScalar **)&Y->x));
215d71ae5a4SJacob Faibussowitsch     break;
216d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_WRITE:
217d71ae5a4SJacob Faibussowitsch     PetscCall(VecRestoreArrayWrite(X, (PetscScalar **)&Y->x));
218d71ae5a4SJacob Faibussowitsch     break;
219d71ae5a4SJacob Faibussowitsch   default:
220d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
221d71ae5a4SJacob Faibussowitsch     break;
222218db3c1SStefano Zampini   }
2233ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
224218db3c1SStefano Zampini }
225218db3c1SStefano Zampini 
MatDenseWrapCholmod(Mat X,PetscInt rw,cholmod_dense * Y)226d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDenseWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y)
227d71ae5a4SJacob Faibussowitsch {
228218db3c1SStefano Zampini   PetscScalar *x;
229218db3c1SStefano Zampini   PetscInt     m, n, lda;
230218db3c1SStefano Zampini 
231218db3c1SStefano Zampini   PetscFunctionBegin;
2329566063dSJacob Faibussowitsch   PetscCall(PetscMemzero(Y, sizeof(*Y)));
233218db3c1SStefano Zampini   switch (rw) {
234d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_READ:
235d71ae5a4SJacob Faibussowitsch     PetscCall(MatDenseGetArrayRead(X, (const PetscScalar **)&x));
236d71ae5a4SJacob Faibussowitsch     break;
237d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_WRITE:
238d71ae5a4SJacob Faibussowitsch     PetscCall(MatDenseGetArrayWrite(X, &x));
239d71ae5a4SJacob Faibussowitsch     break;
240d71ae5a4SJacob Faibussowitsch   default:
241d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
242d71ae5a4SJacob Faibussowitsch     break;
243218db3c1SStefano Zampini   }
2449566063dSJacob Faibussowitsch   PetscCall(MatDenseGetLDA(X, &lda));
2459566063dSJacob Faibussowitsch   PetscCall(MatGetLocalSize(X, &m, &n));
246218db3c1SStefano Zampini 
247218db3c1SStefano Zampini   Y->x     = x;
248218db3c1SStefano Zampini   Y->nrow  = m;
249218db3c1SStefano Zampini   Y->ncol  = n;
250218db3c1SStefano Zampini   Y->nzmax = lda * n;
251218db3c1SStefano Zampini   Y->d     = lda;
252218db3c1SStefano Zampini   Y->xtype = CHOLMOD_SCALAR_TYPE;
253218db3c1SStefano Zampini   Y->dtype = CHOLMOD_DOUBLE;
2543ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
255218db3c1SStefano Zampini }
256218db3c1SStefano Zampini 
MatDenseUnWrapCholmod(Mat X,PetscInt rw,cholmod_dense * Y)257d71ae5a4SJacob Faibussowitsch PetscErrorCode MatDenseUnWrapCholmod(Mat X, PetscInt rw, cholmod_dense *Y)
258d71ae5a4SJacob Faibussowitsch {
259218db3c1SStefano Zampini   PetscFunctionBegin;
260218db3c1SStefano Zampini   switch (rw) {
261d71ae5a4SJacob Faibussowitsch   case GET_ARRAY_READ:
262d71ae5a4SJacob Faibussowitsch     PetscCall(MatDenseRestoreArrayRead(X, (const PetscScalar **)&Y->x));
263d71ae5a4SJacob Faibussowitsch     break;
264218db3c1SStefano Zampini   case GET_ARRAY_WRITE:
265218db3c1SStefano Zampini     /* we don't have MatDenseRestoreArrayWrite */
2669566063dSJacob Faibussowitsch     PetscCall(MatDenseRestoreArray(X, (PetscScalar **)&Y->x));
267218db3c1SStefano Zampini     break;
268d71ae5a4SJacob Faibussowitsch   default:
269d71ae5a4SJacob Faibussowitsch     SETERRQ(PETSC_COMM_SELF, PETSC_ERR_SUP, "Case %" PetscInt_FMT " not handled", rw);
270d71ae5a4SJacob Faibussowitsch     break;
271218db3c1SStefano Zampini   }
2723ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
273d9ca1df4SBarry Smith }
274d9ca1df4SBarry Smith 
MatDestroy_CHOLMOD(Mat F)275d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatDestroy_CHOLMOD(Mat F)
276d71ae5a4SJacob Faibussowitsch {
2776b8f6f9dSBarry Smith   Mat_CHOLMOD *chol = (Mat_CHOLMOD *)F->data;
278641875f9SMatthew G Knepley 
279641875f9SMatthew G Knepley   PetscFunctionBegin;
2803ba16761SJacob Faibussowitsch   if (chol->spqrfact) PetscCallExternal(!SuiteSparseQR_C_free, &chol->spqrfact, chol->common);
2813ba16761SJacob Faibussowitsch   if (chol->factor) PetscCallExternal(!cholmod_X_free_factor, &chol->factor, chol->common);
282a2fc1e05SToby Isaac   if (chol->common->itype == CHOLMOD_INT) {
2833ba16761SJacob Faibussowitsch     PetscCallExternal(!cholmod_finish, chol->common);
284a2fc1e05SToby Isaac   } else {
2853ba16761SJacob Faibussowitsch     PetscCallExternal(!cholmod_l_finish, chol->common);
286a2fc1e05SToby Isaac   }
2879566063dSJacob Faibussowitsch   PetscCall(PetscFree(chol->common));
2889566063dSJacob Faibussowitsch   PetscCall(PetscFree(chol->matrix));
2899566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatFactorGetSolverType_C", NULL));
2902e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorSymbolic_C", NULL));
2912e956fe4SStefano Zampini   PetscCall(PetscObjectComposeFunction((PetscObject)F, "MatQRFactorNumeric_C", NULL));
2929566063dSJacob Faibussowitsch   PetscCall(PetscFree(F->data));
2933ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
294641875f9SMatthew G Knepley }
295641875f9SMatthew G Knepley 
296641875f9SMatthew G Knepley static PetscErrorCode MatSolve_CHOLMOD(Mat, Vec, Vec);
297218db3c1SStefano Zampini static PetscErrorCode MatMatSolve_CHOLMOD(Mat, Mat, Mat);
298641875f9SMatthew G Knepley 
299fcd503bcSBarry Smith /*static const char *const CholmodOrderingMethods[] = {"User","AMD","METIS","NESDIS(default)","Natural","NESDIS(small=20000)","NESDIS(small=4,no constrained)","NESDIS()"};*/
300641875f9SMatthew G Knepley 
MatView_Info_CHOLMOD(Mat F,PetscViewer viewer)301d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatView_Info_CHOLMOD(Mat F, PetscViewer viewer)
302d71ae5a4SJacob Faibussowitsch {
3036b8f6f9dSBarry Smith   Mat_CHOLMOD          *chol = (Mat_CHOLMOD *)F->data;
304641875f9SMatthew G Knepley   const cholmod_common *c    = chol->common;
305641875f9SMatthew G Knepley   PetscInt              i;
306641875f9SMatthew G Knepley 
307641875f9SMatthew G Knepley   PetscFunctionBegin;
3083ba16761SJacob Faibussowitsch   if (F->ops->solve != MatSolve_CHOLMOD) PetscFunctionReturn(PETSC_SUCCESS);
3099566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "CHOLMOD run parameters:\n"));
3109566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPushTab(viewer));
3119566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Pack factors after symbolic factorization: %s\n", chol->pack ? "TRUE" : "FALSE"));
3129566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.dbound            %g  (Smallest absolute value of diagonal entries of D)\n", c->dbound));
3139566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow0             %g\n", c->grow0));
3149566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow1             %g\n", c->grow1));
3159566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.grow2             %u\n", (unsigned)c->grow2));
3169566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.maxrank           %u\n", (unsigned)c->maxrank));
3179566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal_switch %g\n", c->supernodal_switch));
3189566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.supernodal        %d\n", c->supernodal));
3199566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_asis        %d\n", c->final_asis));
3209566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_super       %d\n", c->final_super));
3219566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_ll          %d\n", c->final_ll));
3229566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_pack        %d\n", c->final_pack));
3239566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_monotonic   %d\n", c->final_monotonic));
3249566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.final_resymbol    %d\n", c->final_resymbol));
3259566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.zrelax            [%g,%g,%g]\n", c->zrelax[0], c->zrelax[1], c->zrelax[2]));
3269566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrelax            [%u,%u,%u]\n", (unsigned)c->nrelax[0], (unsigned)c->nrelax[1], (unsigned)c->nrelax[2]));
3279566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.prefer_upper      %d\n", c->prefer_upper));
3289566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.print             %d\n", c->print));
329641875f9SMatthew G Knepley   for (i = 0; i < c->nmethods; i++) {
3309566063dSJacob Faibussowitsch     PetscCall(PetscViewerASCIIPrintf(viewer, "Ordering method %" PetscInt_FMT "%s:\n", i, i == c->selected ? " [SELECTED]" : ""));
3319371c9d4SSatish 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));
332641875f9SMatthew G Knepley   }
3339566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.postorder         %d\n", c->postorder));
3349566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.default_nesdis    %d (use NESDIS instead of METIS for nested dissection)\n", c->default_nesdis));
335641875f9SMatthew G Knepley   /* Statistics */
3369566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.fl                %g (flop count from most recent analysis)\n", c->fl));
3379566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.lnz               %g (fundamental nz in L)\n", c->lnz));
3389566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.anz               %g\n", c->anz));
3399566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.modfl             %g (flop count from most recent update)\n", c->modfl));
3409566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.malloc_count      %g (number of live objects)\n", (double)c->malloc_count));
3419566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_usage      %g (peak memory usage in bytes)\n", (double)c->memory_usage));
3429566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.memory_inuse      %g (current memory usage in bytes)\n", (double)c->memory_inuse));
3439566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_col      %g (number of column reallocations)\n", c->nrealloc_col));
3449566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.nrealloc_factor   %g (number of factor reallocations due to column reallocations)\n", c->nrealloc_factor));
3459566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.ndbounds_hit      %g (number of times diagonal was modified by dbound)\n", c->ndbounds_hit));
3469566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.rowfacfl          %g (number of flops in last call to cholmod_rowfac)\n", c->rowfacfl));
3479566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.aatfl             %g (number of flops to compute A(:,f)*A(:,f)')\n", c->aatfl));
348b9eaa5e8SBarry Smith #if defined(PETSC_USE_SUITESPARSE_GPU)
3499566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPrintf(viewer, "Common.useGPU            %d\n", c->useGPU));
350b9eaa5e8SBarry Smith #endif
3519566063dSJacob Faibussowitsch   PetscCall(PetscViewerASCIIPopTab(viewer));
3523ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
353641875f9SMatthew G Knepley }
354641875f9SMatthew G Knepley 
MatView_CHOLMOD(Mat F,PetscViewer viewer)355d71ae5a4SJacob Faibussowitsch PETSC_INTERN PetscErrorCode MatView_CHOLMOD(Mat F, PetscViewer viewer)
356d71ae5a4SJacob Faibussowitsch {
3579f196a02SMartin Diehl   PetscBool         isascii;
358641875f9SMatthew G Knepley   PetscViewerFormat format;
359641875f9SMatthew G Knepley 
360641875f9SMatthew G Knepley   PetscFunctionBegin;
3619f196a02SMartin Diehl   PetscCall(PetscObjectTypeCompare((PetscObject)viewer, PETSCVIEWERASCII, &isascii));
3629f196a02SMartin Diehl   if (isascii) {
3639566063dSJacob Faibussowitsch     PetscCall(PetscViewerGetFormat(viewer, &format));
36448a46eb9SPierre Jolivet     if (format == PETSC_VIEWER_ASCII_INFO) PetscCall(MatView_Info_CHOLMOD(F, viewer));
365641875f9SMatthew G Knepley   }
3663ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
367641875f9SMatthew G Knepley }
368641875f9SMatthew G Knepley 
MatSolve_CHOLMOD(Mat F,Vec B,Vec X)369d71ae5a4SJacob Faibussowitsch static PetscErrorCode MatSolve_CHOLMOD(Mat F, Vec B, Vec X)
370d71ae5a4SJacob Faibussowitsch {
3716b8f6f9dSBarry Smith   Mat_CHOLMOD  *chol = (Mat_CHOLMOD *)F->data;
372218db3c1SStefano Zampini   cholmod_dense cholB, cholX, *X_handle, *Y_handle = NULL, *E_handle = NULL;
373641875f9SMatthew G Knepley 
374641875f9SMatthew G Knepley   PetscFunctionBegin;
37572f3cec6SPablo Brubeck   if (!F->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
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 
MatMatSolve_CHOLMOD(Mat F,Mat B,Mat X)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 
MatCholeskyFactorNumeric_CHOLMOD(Mat F,Mat A,const MatFactorInfo * info)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;
41672f3cec6SPablo Brubeck   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
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 
MatCholeskyFactorSymbolic_CHOLMOD(Mat F,Mat A,IS perm,const MatFactorInfo * info)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;
44772f3cec6SPablo Brubeck   F->ops->choleskyfactornumeric = MatCholeskyFactorNumeric_CHOLMOD;
44872f3cec6SPablo Brubeck   if (!A->rmap->n) PetscFunctionReturn(PETSC_SUCCESS);
44926cc229bSBarry Smith   /* Set options to F */
45026cc229bSBarry Smith   PetscCall(CholmodSetOptions(F));
45126cc229bSBarry Smith 
4529566063dSJacob Faibussowitsch   PetscCall((*chol->Wrap)(A, PETSC_FALSE, &cholA, &aijalloc, &valloc));
453641875f9SMatthew G Knepley   static_F = F;
454641875f9SMatthew G Knepley   if (chol->factor) {
455d0609cedSBarry Smith     err = !cholmod_X_resymbol(&cholA, fset, fsize, (int)chol->pack, chol->factor, chol->common);
456d0609cedSBarry Smith     PetscCheck(!err, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed with status %d", chol->common->status);
457641875f9SMatthew G Knepley   } else if (perm) {
458641875f9SMatthew G Knepley     const PetscInt *ip;
4599566063dSJacob Faibussowitsch     PetscCall(ISGetIndices(perm, &ip));
460641875f9SMatthew G Knepley     chol->factor = cholmod_X_analyze_p(&cholA, (PetscInt *)ip, fset, fsize, chol->common);
46128b400f6SJacob Faibussowitsch     PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using PETSc ordering with status %d", chol->common->status);
4629566063dSJacob Faibussowitsch     PetscCall(ISRestoreIndices(perm, &ip));
463641875f9SMatthew G Knepley   } else {
464641875f9SMatthew G Knepley     chol->factor = cholmod_X_analyze(&cholA, chol->common);
46528b400f6SJacob Faibussowitsch     PetscCheck(chol->factor, PetscObjectComm((PetscObject)F), PETSC_ERR_LIB, "CHOLMOD analysis failed using internal ordering with status %d", chol->common->status);
466641875f9SMatthew G Knepley   }
467641875f9SMatthew G Knepley 
4689566063dSJacob Faibussowitsch   if (aijalloc) PetscCall(PetscFree2(cholA.p, cholA.i));
4699566063dSJacob Faibussowitsch   if (valloc) PetscCall(PetscFree(cholA.x));
4703ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
471641875f9SMatthew G Knepley }
472641875f9SMatthew G Knepley 
MatFactorGetSolverType_seqsbaij_cholmod(Mat A,MatSolverType * type)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 
MatGetInfo_CHOLMOD(Mat F,MatInfoType flag,MatInfo * info)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:
5301d27aa22SBarry Smith    CHOLMOD is part of SuiteSparse <http://faculty.cse.tamu.edu/davis/suitesparse.html>
531a364b7d2SBarry Smith 
5321cc06b55SBarry Smith .seealso: [](ch_matrices), `Mat`, `PCCHOLESKY`, `PCFactorSetMatSolverType()`, `MatSolverType`
533641875f9SMatthew G Knepley M*/
534b2573a8aSBarry Smith 
MatGetFactor_seqsbaij_cholmod(Mat A,MatFactorType ftype,Mat * F)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));
54303e5aca4SStefano Zampini   *F = NULL;
54403e5aca4SStefano Zampini   if (bs != 1) {
54503e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "CHOLMOD only supports block size=1.\n"));
54603e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
54703e5aca4SStefano Zampini   }
548*fc2fb351SPierre Jolivet   if (PetscDefined(USE_COMPLEX) && A->hermitian != PETSC_BOOL3_TRUE) {
54903e5aca4SStefano Zampini     PetscCall(PetscInfo(A, "Only for Hermitian matrices.\n"));
55003e5aca4SStefano Zampini     PetscFunctionReturn(PETSC_SUCCESS);
55103e5aca4SStefano Zampini   }
552641875f9SMatthew G Knepley   /* Create the factorization matrix F */
5539566063dSJacob Faibussowitsch   PetscCall(MatCreate(PetscObjectComm((PetscObject)A), &B));
5549566063dSJacob Faibussowitsch   PetscCall(MatSetSizes(B, PETSC_DECIDE, PETSC_DECIDE, m, n));
5559566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy("cholmod", &((PetscObject)B)->type_name));
5569566063dSJacob Faibussowitsch   PetscCall(MatSetUp(B));
5574dfa11a4SJacob Faibussowitsch   PetscCall(PetscNew(&chol));
55826fbe8dcSKarl Rupp 
559641875f9SMatthew G Knepley   chol->Wrap = MatWrapCholmod_seqsbaij;
5606b8f6f9dSBarry Smith   B->data    = chol;
561641875f9SMatthew G Knepley 
562218db3c1SStefano Zampini   B->ops->getinfo                = MatGetInfo_CHOLMOD;
563641875f9SMatthew G Knepley   B->ops->view                   = MatView_CHOLMOD;
564641875f9SMatthew G Knepley   B->ops->choleskyfactorsymbolic = MatCholeskyFactorSymbolic_CHOLMOD;
565641875f9SMatthew G Knepley   B->ops->destroy                = MatDestroy_CHOLMOD;
5669566063dSJacob Faibussowitsch   PetscCall(PetscObjectComposeFunction((PetscObject)B, "MatFactorGetSolverType_C", MatFactorGetSolverType_seqsbaij_cholmod));
567641875f9SMatthew G Knepley   B->factortype   = MAT_FACTOR_CHOLESKY;
568218db3c1SStefano Zampini   B->assembled    = PETSC_TRUE;
569641875f9SMatthew G Knepley   B->preallocated = PETSC_TRUE;
570641875f9SMatthew G Knepley 
5719566063dSJacob Faibussowitsch   PetscCall(CholmodStart(B));
57200c67f3bSHong Zhang 
5739566063dSJacob Faibussowitsch   PetscCall(PetscFree(B->solvertype));
5749566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATSOLVERCHOLMOD, &B->solvertype));
575f73b0415SBarry Smith   B->canuseordering = PETSC_TRUE;
5769566063dSJacob Faibussowitsch   PetscCall(PetscStrallocpy(MATORDERINGEXTERNAL, (char **)&B->preferredordering[MAT_FACTOR_CHOLESKY]));
577641875f9SMatthew G Knepley   *F = B;
5783ba16761SJacob Faibussowitsch   PetscFunctionReturn(PETSC_SUCCESS);
579641875f9SMatthew G Knepley }
580